API Reference

The table below lists the main entry points exported by Geodynamo. The listing is automatically generated during the documentation build.

Main Module

Geodynamo.ERK2CacheType
ERK2Cache{T}

Cached data structure for Exponential 2nd Order Runge-Kutta method. Stores precomputed matrix exponentials and φ functions for each spherical harmonic mode.

source
Geodynamo.FieldCombinerType
FieldCombiner{T}

Structure for combining distributed Geodynamo.jl output files using the consistent field structures and parameter system.

source
Geodynamo._get_influence_cacheMethod
_get_influence_cache(domain::RadialDomain, dr_matrix)

Get or create cached influence functions and matrix for given domain. Note: drmatrix should be a BandedMatrix (defined in linearalgebra.jl)

source
Geodynamo.add_internal_sources_local!Method
add_internal_sources_local!(field::AbstractScalarField{T}, domain::RadialDomain) where T

Add volumetric sources (completely local operation). This works for any scalar field with radial source profile.

source
Geodynamo.apply_flux_bc_direct!Method
apply_flux_bc_direct!(spec_real, spec_imag, local_lm, lm_idx,
                     field::AbstractScalarField, domain, r_range)

Direct modification of boundary values to approximately satisfy flux BC. This is the simplest but least accurate method.

source
Geodynamo.apply_flux_bc_influence_matrix!Method
apply_flux_bc_influence_matrix!(spec_real, spec_imag, local_lm, lm_idx,
                               apply_inner, apply_outer, field::AbstractScalarField,
                               domain, r_range)

Apply flux boundary conditions using the influence matrix method.

source
Geodynamo.apply_flux_bc_spectral!Method
apply_flux_bc_spectral_complete!(temp_field, domain)

Complete implementation of flux boundary conditions in spectral space. This modifies the spectral coefficients to satisfy ∂T/∂r = prescribed_flux.

source
Geodynamo.apply_flux_bc_tau!Method
apply_flux_bc_tau!(spec_real, spec_imag, local_lm, lm_idx, 
                   apply_inner, apply_outer, field::AbstractScalarField, 
                   domain, r_range)

Apply flux boundary conditions using the tau method. This is the generalized version that works with any scalar field.

source
Geodynamo.apply_geometric_factors_spectral!Method
apply_geometric_factors_spectral!(field::AbstractScalarField{T}, domain::RadialDomain) where T

Apply geometric factors (1/r, 1/(r sin θ)) in spectral space. For gradients in spherical coordinates. This is generic.

source
Geodynamo.apply_magnetic_boundary_conditions!Method
apply_magnetic_boundary_conditions!(fields; time_index=nothing)

Refresh magnetic boundary data from the BoundaryConditions subsystem and enforce the corresponding Dirichlet values in spectral space.

source
Geodynamo.apply_scalar_flux_bc_spectral!Method
apply_scalar_flux_bc_spectral!(field::AbstractScalarField{T}, domain::RadialDomain;
                               method::Symbol=:tau) where T

Apply flux boundary conditions to a scalar field in spectral space. Methods available: :tau (most robust), :influence_matrix, :direct (simplest).

source
Geodynamo.apply_velocity_boundary_conditions!Method
apply_velocity_boundary_conditions!(fields; time_index=nothing)

Refresh velocity boundary data from the BoundaryConditions subsystem and immediately enforce Dirichlet constraints in spectral space.

source
Geodynamo.batch_convert_directoryMethod
batch_convert_directory(input_dir::String, output_dir::String = "";
                       pattern::String = "combined_global_time_",
                       precision::Type{T} = Float64,
                       compute_diagnostics::Bool = true) where T

Convert all spectral files matching a pattern in a directory.

source
Geodynamo.batch_shtnskit_transforms!Method
batch_shtnskit_transforms!(specs::Vector{SHTnsSpectralField{T}},
                          physs::Vector{SHTnsPhysicalField{T}}) where T

Batch process multiple transforms using SHTnsKit with PencilArrays.

source
Geodynamo.batch_spectral_to_physical!Method
batch_spectral_to_physical!(specs, physs)

Compatibility wrapper that calls batch_shtnskit_transforms! for batched spectral→physical transforms using SHTnsKit with PencilArrays/MPI.

source
Geodynamo.batch_write_spectral_fields!Method
batch_write_spectral_fields!(nc_file, fields::Dict, config::OutputConfig, field_info::FieldInfo)

Write multiple spectral fields in batched operations for better I/O performance

source
Geodynamo.build_banded_AMethod
build_banded_A(T, domain, diffusivity, l) -> BandedMatrix{T}

Construct banded A = ν*(d²/dr² + (2/r)d/dr − l(l+1)/r²) in banded storage.

source
Geodynamo.build_influence_matrixMethod
build_influence_matrix(G_inner, G_outer, dr_matrix, domain)

Construct a 2×2 matrix mapping influence amplitudes to boundary flux errors. Rows correspond to (inner, outer) boundaries; columns to (inner, outer) influence functions.

source
Geodynamo.build_rhs_cnab2!Method
build_rhs_cnab2!(rhs, un, nl, nl_prev, dt, matrices)

Build RHS for CNAB2 IMEX: rhs = un/dt + (1-θ)·L·un + (3/2)·nl − (1/2)·nl_prev, where θ = matrices.theta and L is the diffusivity-scaled linear operator.

source
Geodynamo.build_theta_derivative_matrixMethod
build_theta_derivative_matrix(::Type{T}, config::SHTnsKitConfig) where T

Build sparse matrix for θ-derivatives in spectral space. This matrix couples different l modes with the same m.

source
Geodynamo.combine_distributed_timeMethod
combine_distributed_time(output_dir::String, time::Float64; 
                        config::CombinerConfig = create_combiner_config(),
                        precision::Type{T} = Float64) where T

Main function to combine distributed files for a specific time.

source
Geodynamo.combine_time_seriesMethod
combine_time_series(output_dir::String, time_range::Tuple{Float64,Float64}; 
                   config::CombinerConfig = create_combiner_config(),
                   precision::Type{T} = Float64) where T

Combine multiple time points into a time series.

source
Geodynamo.compute_all_gradients_spectral!Method
compute_all_gradients_spectral!(field::AbstractScalarField{T}, domain::RadialDomain) where T

Compute all gradient components (θ, φ, r) in spectral space. This is the main driver function that works for any scalar field.

source
Geodynamo.compute_boundary_fluxesMethod
compute_boundary_fluxes(profile::Vector{T}, dr_matrix, domain::RadialDomain) where T

Compute flux (dT/dr) at both boundaries using the derivative matrix. Note: drmatrix should be a BandedMatrix{T} (defined in linearalgebra.jl)

source
Geodynamo.compute_nonlinear!Method
compute_nonlinear!(parallelizer, temp_field, vel_fields, domain)

Advanced CPU computation with SIMD, NUMA, and task-based parallelism.

source
Geodynamo.compute_phi1_functionMethod
compute_phi1_function(A, expA)

Compute φ1(A) = (exp(A) - I) / A efficiently with comprehensive error handling. Uses series expansion for small ||A|| to avoid numerical issues.

source
Geodynamo.compute_phi2_functionMethod
compute_phi2_function(A, expA)

Compute φ2(A) = (exp(A) - I - A) / A² efficiently with comprehensive error handling. Uses series expansion for small ||A|| to avoid numerical issues.

source
Geodynamo.compute_phi_gradient_spectral!Method
compute_phi_gradient_spectral!(field::AbstractScalarField{T}) where T

Compute ∂field/∂φ using spherical harmonic properties (local operation). This is generic and works for any scalar field.

source
Geodynamo.compute_radial_gradient_spectral!Method
compute_radial_gradient_spectral!(field::AbstractScalarField{T}, domain::RadialDomain) where T

Compute ∂field/∂r using banded matrix derivative operator in spectral space (local operation). This uses the pre-computed derivative matrices from the field for optimal accuracy and efficiency.

source
Geodynamo.compute_scalar_advection_local!Method
compute_scalar_advection_local!(field::AbstractScalarField{T}, vel_fields) where T

Compute -u·∇field in physical space (completely local operation). This works for any scalar field.

source
Geodynamo.compute_tau_coefficients_bothMethod
compute_tau_coefficients_both(flux_error_inner::T, flux_error_outer::T, domain::RadialDomain) where T

Compute tau polynomial coefficients for both boundaries. Uses highest two Chebyshev modes as tau functions.

source
Geodynamo.compute_tau_coefficients_innerMethod
compute_tau_coefficients_inner(flux_error::T, domain::RadialDomain) where T

Compute tau coefficient for inner boundary only. Uses a single tau function that doesn't affect outer boundary.

source
Geodynamo.compute_tau_coefficients_outerMethod
compute_tau_coefficients_outer(flux_error::T, domain::RadialDomain) where T

Compute tau coefficient for outer boundary only. Uses a single tau function that doesn't affect inner boundary.

source
Geodynamo.compute_theta_gradient_spectral!Method
compute_theta_gradient_spectral!(field::AbstractScalarField{T}) where T

Compute ∂field/∂θ using spherical harmonic recurrence relations (local operation). This is generic and works for any scalar field.

source
Geodynamo.compute_theta_recurrence_coefficientsMethod
compute_theta_recurrence_coefficients(::Type{T}, config::SHTnsKitConfig) where T

Pre-compute recurrence coefficients for θ-derivatives. Store as [nlm, 2] matrix: [:, 1] for l-1 coupling, [:, 2] for l+1 coupling

source
Geodynamo.convert_spectral_fileMethod
convert_spectral_file(input_filename::String, output_filename::String = "";
                     precision::Type{T} = Float64, 
                     compute_diagnostics::Bool = true) where T

Main function to convert a single spectral data file to physical space.

source
Geodynamo.convert_to_physical!Method
convert_to_physical!(converter::SpectralToPhysicalConverter{T}) where T

Convert all loaded spectral data to physical space using SHTns transforms.

source
Geodynamo.create_enhanced_field_variables!Method
create_enhanced_field_variables!(nc_file, field_info::FieldInfo, config::OutputConfig, 
                                available_fields::Vector{String})

Enhanced field variable setup that leverages SHTns configuration

source
Geodynamo.create_erk2_cacheMethod
create_erk2_cache(config, domain, diffusivity, dt; use_krylov=false, m=20, tol=1e-8)

Create ERK2 cache with precomputed matrix functions for all spherical harmonic modes.

source
Geodynamo.create_etd_cacheMethod
create_etd_cache(config, domain, diffusivity, dt) -> ETDCache

Build per-l exponential cache for the linear operator Al = diffusivity*(d²/dr² + (2/r)d/dr − l(l+1)/r²). Computes exp(dt Al) and phi1(dt A_l) via dense methods. Single-rank recommended.

source
Geodynamo.create_field_combinerMethod
create_field_combiner(output_dir::String, time::Float64; 
                     precision::Type{T} = Float64) where T

Create a field combiner by reading configuration from distributed files.

source
Geodynamo.create_pencil_topologyMethod
create_pencil_topology(shtns_config; optimize=true)

Create enhanced pencil decomposition for SHTns grids. Supports both 1D and 2D decompositions with automatic optimization. Accepts an object with fields nlat, nlon, and nlm (e.g., SHTnsKitConfig).

source
Geodynamo.create_shtns_aware_output_configMethod
create_shtns_aware_output_config(shtns_config::SHTnsKitConfig, pencils::NamedTuple; 
                                base_config::OutputConfig = default_config())

Create output configuration that integrates with SHTns configuration and pencil decomposition

source
Geodynamo.create_shtnskit_configMethod
create_shtnskit_config(; lmax::Int, mmax::Int=lmax, nlat::Int=lmax+2,
                       nlon::Int=max(2*lmax+1, 4), optimize_decomp::Bool=true) -> SHTnsKitConfig

Create SHTnsKit configuration with MPI parallelization using PencilArrays. This creates proper integration with ../SHTnsKit.jl, PencilArrays, and PencilFFTs.

source
Geodynamo.create_spectral_converterMethod
create_spectral_converter(filename::String; precision::Type{T} = Float64) where T

Create a spectral to physical converter by reading configuration from a NetCDF file.

source
Geodynamo.eab2_update!Method
eab2_update!(u, nl, nl_prev, etd, config)

Apply EAB2 update per (l,m): u^{n+1} = E u^n + dtphi1(3/2 nl^n − 1/2 nl^{n−1}).

source
Geodynamo.eab2_update_krylov!Method
eab2_update_krylov!(u, nl, nl_prev, domain, diffusivity, config, dt; m=20, tol=1e-8)

EAB2 update using Krylov exp/φ1 actions and banded LU for φ1.

source
Geodynamo.eab2_update_krylov_cached!Method
eab2_update_krylov_cached!(u, nl, nl_prev, alu_map, domain, ν, config, dt; m=20, tol=1e-8)

Same as eab2updatekrylov!, but reuses cached banded A and LU per l.

source
Geodynamo.erk2_stage_residual_statsMethod
erk2_stage_residual_stats(buffers) -> NamedTuple

Compute diagnostic statistics for the difference between stage nonlinear terms and the base-step nonlinear terms.

source
Geodynamo.exchange_dimension_halos!Method
exchange_dimension_halos!(data::Array, pencil::Pencil, dim::Int, 
                          halo_width::Int, boundaries::Symbol, comm::MPI.Comm)

Perform halo exchange along a specific dimension using MPI point-to-point communication.

source
Geodynamo.exp_action_krylovMethod
exp_action_krylov(Aop!, v, dt; m=20, tol=1e-8) -> y ≈ exp(dt A) v

Simple Arnoldi-based approximation of the exponential action.

source
Geodynamo.extract_local_radial_profile!Method
extract_local_radial_profile!(profile, data, local_lm, nr, r_range)

In-place version to avoid allocations; writes the local radial line into profile for the given local_lm using the provided r_range.

source
Geodynamo.find_distributed_filesFunction
find_distributed_files(output_dir::String, time::Float64, filename_prefix::String = "geodynamo")

Find all distributed files for a specific simulation time.

source
Geodynamo.get_commMethod
get_comm()

Get MPI communicator, initializing MPI if needed. Provides thread-safe lazy initialization.

source
Geodynamo.get_dimension_neighborsMethod
get_dimension_neighbors(pencil::Pencil, dim::Int, boundaries::Symbol) -> (Int, Int)

Get left and right neighbor process ranks for a given dimension. Returns (leftneighbor, rightneighbor) where MPI.MPIPROCNULL indicates no neighbor.

source
Geodynamo.get_eab2_alu_cache!Method
get_eab2_alu_cache!(caches, key, ν, T, domain) -> Dict{Int,Tuple{BandedMatrix{T},BandedLU{T}}}

Retrieve or initialize a cache mapping l -> (Abanded, LU(Abanded)) for EAB2. Reinitializes if ν or nr changed.

source
Geodynamo.get_erk2_cache!Method
get_erk2_cache!(caches, key, diffusivity, config, domain, dt; use_krylov=false)

Retrieve or create ERK2 cache with automatic invalidation when parameters change.

source
Geodynamo.get_flux_valueMethod
get_flux_value(lm_idx::Int, boundary::Int, field::AbstractScalarField)

Get prescribed flux value for given mode and boundary from scalar field. This is a generic version that works with any scalar field.

source
Geodynamo.infer_pencils_from_transpose_nameMethod
infer_pencils_from_transpose_name(name::Symbol, plan) -> (src_pencil, dest_pencil)

Infer source and destination pencils from transpose operation name and plan. For TransposeOperator, we need to extract pencil info from context since the operator itself doesn't expose src/dest pencils directly.

source
Geodynamo.list_available_timesFunction
list_available_times(output_dir::String, filename_prefix::String = "geodynamo")

List all available simulation times in the output directory.

source
Geodynamo.load_parametersFunction
load_parameters(config_file::String = "")

Load parameters from a configuration file. If no file is specified, loads from the default config/default_params.jl file.

Arguments

  • config_file::String: Path to parameter file (optional)

Returns

  • GeodynamoParameters: Loaded parameters
source
Geodynamo.load_single_file!Method
load_single_file!(combiner::FieldCombiner{T}, filename::String) where T

Load data from a single file (no combination needed).

source
Geodynamo.load_spectral_coefficients!Method
load_spectral_coefficients!(field::SHTnsSpectralField{T}, nc_file, 
                            real_var_name::String, imag_var_name::String,
                            l_values::Vector{Int}, m_values::Vector{Int}, 
                            r_spectral::Vector{Float64}) where T

Load spectral coefficients from NetCDF file into a spectral field structure.

source
Geodynamo.load_spectral_data!Method
load_spectral_data!(converter::SpectralToPhysicalConverter{T}, filename::String) where T

Load spectral field data from NetCDF file into Geodynamo field structures.

source
Geodynamo.main_combine_timeMethod
main_combine_time(output_dir::String, time::Float64; kwargs...)

Main entry point for combining distributed files for a single time.

source
Geodynamo.map_spectral_coefficients!Method
map_spectral_coefficients!(global_tor_real, global_tor_imag, global_pol_real, global_pol_imag,
                          local_tor_real, local_tor_imag, local_pol_real, local_pol_imag,
                          local_l, local_m, lm_mapping)

Map local spectral coefficients to global arrays using (l,m) index mapping.

source
Geodynamo.modify_for_flux_inner!Method
modify_for_flux_inner!(spec_real, spec_imag, local_lm, prescribed_flux,
                      dr_matrix, domain, r_range)

Modify coefficients near inner boundary to approximate flux condition. Uses low-order extrapolation.

source
Geodynamo.modify_for_flux_outer!Method
modify_for_flux_outer!(spec_real, spec_imag, local_lm, prescribed_flux,
                      dr_matrix, domain, r_range)

Modify coefficients near outer boundary to approximate flux condition. Uses low-order extrapolation.

source
Geodynamo.optimize_erk2_transforms!Method
optimize_erk2_transforms!(config::SHTnsKitConfig)

Optimize SHTnsKit transforms for ERK2 timestepping with PencilFFTs. This function pre-warms transform plans and optimizes memory layout.

source
Geodynamo.optimize_process_topologyMethod
optimize_process_topology(nprocs::Int, dims::Tuple{Int,Int,Int})

Find optimal 2D process grid for given number of processes and problem dimensions. Minimizes communication volume.

source
Geodynamo.phi1_action_krylovMethod
phi1_action_krylov(BA, LU_A, v, dt; m=20, tol=1e-8) -> y ≈ φ1(dt A) v

Compute φ1(dt A) v = A^{-1}[(exp(dt A) − I) v]/dt using Krylov exp(action) and banded solve.

source
Geodynamo.print_pencil_axesMethod
print_pencil_axes(pencils)

Print the axes_local tuple for each pencil, showing the local index ranges for all three axes. This helps verify which axes are distributed (those with nontrivial subranges across ranks) and which axis is contiguous locally.

source
Geodynamo.run_simulation!Method
run_master_simulation!(state::SimulationState{T})

Run geodynamo simulation with maximum CPU parallelization optimizations.

source
Geodynamo.save_combined_fieldsMethod
save_combined_fields(combiner::FieldCombiner{T}, output_filename::String; 
                    config::CombinerConfig = create_combiner_config()) where T

Save combined fields to NetCDF file using Geodynamo I/O system.

source
Geodynamo.save_combined_time_seriesMethod
save_combined_time_series(combiners::Vector{FieldCombiner{T}}, output_dir::String, 
                         filename_prefix::String; 
                         config::CombinerConfig = create_combiner_config()) where T

Save combined time series to a single NetCDF file.

source
Geodynamo.save_physical_fieldsMethod
save_physical_fields(converter::SpectralToPhysicalConverter{T}, output_filename::String) where T

Save converted physical space fields to NetCDF file using parallel I/O.

source
Geodynamo.set_velocity_workspace!Method
set_velocity_workspace!(ws)

Register a global VelocityWorkspace to be used by velocity kernels when available. Pass nothing to disable and fall back to internal buffers.

source
Geodynamo.shtnskit_spectral_to_physical!Method
shtnskit_spectral_to_physical!(spec::SHTnsSpectralField{T},
                              phys::SHTnsPhysicalField{T}) where T

Transform from spectral to physical space using SHTnsKit with PencilArrays/PencilFFTs.

source
Geodynamo.shtnskit_vector_analysis!Method
shtnskit_vector_analysis!(vec_phys::SHTnsVectorField{T},
                         tor_spec::SHTnsSpectralField{T}, 
                         pol_spec::SHTnsSpectralField{T}) where T

Vector analysis using SHTnsKit with PencilArrays.

source
Geodynamo.shtnskit_vector_synthesis!Method
shtnskit_vector_synthesis!(tor_spec::SHTnsSpectralField{T}, 
                          pol_spec::SHTnsSpectralField{T},
                          vec_phys::SHTnsVectorField{T}) where T

Vector synthesis using SHTnsKit spheroidal-toroidal decomposition with PencilArrays.

source
Geodynamo.synchronize_halos!Method
synchronize_halos!(arr::PencilArray; halo_width::Int=1, boundaries::Symbol=:all)

Synchronize ghost/halo regions for parallel finite difference computations.

Parameters

  • arr: PencilArray to synchronize
  • halo_width: Width of halo region (number of ghost cells), default 1
  • boundaries: Boundary handling (:all, :periodic, :nonperiodic), default :all

Notes

This function implements MPI point-to-point communication to exchange boundary data between neighboring processes in each decomposed dimension. It's essential for maintaining accuracy in finite difference stencil operations near subdomain boundaries.

For spectral methods, explicit halo exchange may not be needed as the global nature of spectral transforms handles boundary coupling through transpose operations.

source
Geodynamo.test_halo_exchangeMethod
test_halo_exchange(pencil::Pencil, ::Type{T}=Float64; halo_width::Int=1, verbose::Bool=true) where T

Test halo exchange functionality by creating a test array with rank-based values. Returns true if halo exchange is working correctly.

Example

pencils = create_pencil_topology(shtns_config)
test_halo_exchange(pencils.θ, Float64; halo_width=1, verbose=true)
source
Geodynamo.write_grid_file!Method
write_grid_file!(config::OutputConfig, field_info::FieldInfo,
                shtns_config::Union{SHTnsKitConfig,Nothing},
                metadata::Dict{String,Any})

Write a separate grid file containing all coordinate and grid information. This is written only once by rank 0 at the start of the simulation.

source
Geodynamo.write_spectral_data_parallel!Method
write_spectral_data_parallel!(nc_file, real_name, imag_name, real_data, imag_data, field_info)

Write spectral data using parallel I/O strategies based on pencil decomposition

source

Boundary Conditions

Geodynamo.BoundaryConditions.apply_composition_boundaries!Method
apply_composition_boundaries!(comp_field, boundary_set::BoundaryConditionSet; time::Float64=0.0)

Apply composition boundary conditions from a BoundaryConditionSet to a composition field. This is a convenience wrapper that loads the boundary set and applies it.

source
Geodynamo.BoundaryConditions.apply_temperature_boundaries!Method
apply_temperature_boundaries!(temp_field, boundary_set::BoundaryConditionSet; time::Float64=0.0)

Apply temperature boundary conditions from a BoundaryConditionSet to a temperature field. This is a convenience wrapper that loads the boundary set and applies it.

source
Geodynamo.BoundaryConditions.bilinear_interpolateMethod
bilinear_interpolate(data::Matrix{T}, theta_idx::Tuple{Int, Int}, phi_idx::Tuple{Int, Int},
                    theta_weights::Tuple{T, T}, phi_weights::Tuple{T, T}) where T

Perform bilinear interpolation on a 2D data array. Handles periodic boundary conditions in phi direction.

source
Geodynamo.BoundaryConditions.compute_boundary_condition_residualMethod
compute_boundary_condition_residual(field, field_type::FieldType)

Compute residual to check how well boundary conditions are satisfied.

Returns a measure of how much the current field violates the boundary conditions. Useful for monitoring solution quality and debugging.

source
Geodynamo.BoundaryConditions.create_layered_composition_boundaryMethod
create_layered_composition_boundary(config, layer_specs::Vector{Tuple{Real, Real, Real}})

Create layered composition boundary conditions.

Arguments

  • layer_specs: Vector of (colatitudestart, colatitudeend, composition) tuples

Example

# Create three-layer composition structure
layer_specs = [
    (0.0, π/3, 0.8),     # High composition in top layer
    (π/3, 2π/3, 0.4),    # Medium composition in middle layer  
    (2π/3, π, 0.1)       # Low composition in bottom layer
]
source
Geodynamo.BoundaryConditions.create_layered_temperature_boundaryMethod
create_layered_temperature_boundary(config, layer_specs::Vector{Tuple{Real, Real, Real}})

Create layered temperature boundary conditions.

Arguments

  • layer_specs: Vector of (colatitudestart, colatitudeend, temperature) tuples

Example

# Create three-layer temperature structure
layer_specs = [
    (0.0, π/3, 1000.0),    # High temperature in top layer
    (π/3, 2π/3, 500.0),   # Medium temperature in middle layer  
    (2π/3, π, 100.0)      # Low temperature in bottom layer
]
source
Geodynamo.BoundaryConditions.create_programmatic_boundaryFunction
create_programmatic_boundary(pattern::Symbol, config, amplitude::Real=1.0;
                            parameters::Dict=Dict(), description::String="")

Create programmatically generated boundary conditions.

Available patterns:

  • :uniform - Uniform value
  • :y11 - Y₁₁ spherical harmonic
  • :plume - Gaussian plume pattern
  • :hemisphere - Hemispherical pattern
  • :dipole - Dipolar pattern (Y₁₀)
  • :quadrupole - Quadrupolar pattern
  • :custom - User-defined function
source
Geodynamo.BoundaryConditions.create_programmatic_magnetic_boundaryFunction
create_programmatic_magnetic_boundary(pattern::Symbol, config, amplitude::Real=1.0; 
                                    parameters::Dict=Dict())

Create programmatically generated magnetic boundary conditions.

Available patterns:

  • :insulating - Insulating boundary (Br = 0, ∂Btan/∂r = 0)
  • :perfect_conductor - Perfect conductor (B_tan = 0)
  • :dipole - Dipolar magnetic field pattern
  • :quadrupole - Quadrupolar magnetic field pattern
  • :potential_field - Potential field from spherical harmonic coefficients
  • :uniform_field - Uniform magnetic field
  • :custom - User-defined magnetic field function
source
Geodynamo.BoundaryConditions.create_programmatic_velocity_boundaryFunction
create_programmatic_velocity_boundary(pattern::Symbol, config, amplitude::Real=1.0; 
                                    parameters::Dict=Dict())

Create programmatically generated velocity boundary conditions.

Available patterns:

  • :no_slip - Zero velocity at boundary (amplitude ignored)
  • :stress_free - Zero stress at boundary (amplitude ignored)
  • :uniform_rotation - Uniform rotation with angular velocity amplitude
  • :differential_rotation - Differential rotation pattern
  • :zonal_flow - Zonal (east-west) flow pattern
  • :meridional_flow - Meridional (north-south) flow pattern
  • :custom - User-defined velocity function
source
Geodynamo.BoundaryConditions.enforce_boundary_conditions_in_solution!Method
enforce_boundary_conditions_in_solution!(solution, state, field_type::FieldType)

Enforce boundary conditions in the solution vector after timestepping.

This function ensures that the solution satisfies the boundary conditions after each timestep, which may be necessary for certain discretization schemes.

source
Geodynamo.BoundaryConditions.enforce_magnetic_boundary_constraints!Method
enforce_magnetic_boundary_constraints!(magnetic_field, bc_type::Symbol)

Enforce magnetic boundary condition constraints based on physics.

Boundary condition types:

  • :insulating - Insulating boundary: Jnormal = 0, ∇×Btangential = 0 (B_r ≠ 0)
  • :perfect_conductor - Perfect conductor: Btangential = 0, Br free
  • :potential_field - Potential field matching at boundary
source
Geodynamo.BoundaryConditions.enforce_velocity_boundary_constraints!Function
enforce_velocity_boundary_constraints!(velocity_field, bc_type::Symbol=:no_slip)

Enforce specific velocity boundary constraints based on boundary condition type. Uses proper QST (Q-spheroidal-toroidal) decomposition from SHTnsKit.

Arguments

  • velocity_field: Velocity field structure with QST components
  • bc_type: Type of boundary condition (:noslip, :stressfree, :custom)

Boundary Condition Mapping:

  • No-slip: vr = vθ = v_φ = 0 → Q = S = T = 0 at boundaries
  • Stress-free: v_r = 0, tangential stress = 0 → Q = 0, ∂S/∂r = ∂T/∂r = 0
  • Impermeable: v_r = 0 → Q = 0 at boundaries (S, T unconstrained)
source
Geodynamo.BoundaryConditions.find_grid_indicesMethod
find_grid_indices(coords::Vector{T}, target::T; is_periodic::Bool=false) where T

Find the two surrounding indices in a coordinate array for interpolation. Handles periodic coordinates (e.g., longitude) if is_periodic=true.

source
Geodynamo.BoundaryConditions.initialize_boundary_conditions!Method
initialize_boundary_conditions!(field, field_type::FieldType, config)

Initialize boundary condition support for a field structure.

Adds the necessary boundary condition fields to existing field structures without breaking compatibility with existing code.

source
Geodynamo.BoundaryConditions.load_composition_boundary_conditions!Method
load_composition_boundary_conditions!(comp_field, boundary_specs::Dict)

Load composition boundary conditions from various sources.

Arguments

  • comp_field: SHTnsCompositionField structure
  • boundary_specs: Dictionary specifying boundary sources

Examples

# NetCDF files for both boundaries
boundary_specs = Dict(
    :inner => "cmb_composition.nc",
    :outer => "surface_composition.nc"
)

# Mixed NetCDF and programmatic
boundary_specs = Dict(
    :inner => "cmb_composition.nc", 
    :outer => (:uniform, 0.1)  # 10% concentration
)
source
Geodynamo.BoundaryConditions.load_magnetic_boundary_conditions!Method
load_magnetic_boundary_conditions!(magnetic_field, boundary_specs::Dict)

Load magnetic field boundary conditions from various sources.

Arguments

  • magnetic_field: SHTnsMagneticField structure
  • boundary_specs: Dictionary specifying boundary sources

Examples

# Insulating inner, potential field outer
boundary_specs = Dict(
    :inner => (:insulating, 0.0),
    :outer => (:potential_field, "geomagnetic_coefficients.nc")
)

# Perfect conductor boundaries
boundary_specs = Dict(
    :inner => (:perfect_conductor, 0.0),
    :outer => (:perfect_conductor, 0.0)
)

# NetCDF files for both boundaries
boundary_specs = Dict(
    :inner => "cmb_magnetic.nc",
    :outer => "surface_magnetic.nc"
)
source
Geodynamo.BoundaryConditions.load_temperature_boundary_conditions!Method
load_temperature_boundary_conditions!(temp_field, boundary_specs::Dict)

Load temperature boundary conditions from various sources.

Arguments

  • temp_field: SHTnsTemperatureField structure
  • boundary_specs: Dictionary specifying boundary sources

Examples

# NetCDF files for both boundaries
boundary_specs = Dict(
    :inner => "cmb_temperature.nc",
    :outer => "surface_temperature.nc"
)

# Mixed NetCDF and programmatic
boundary_specs = Dict(
    :inner => "cmb_temperature.nc", 
    :outer => (:uniform, 300.0)
)
source
Geodynamo.BoundaryConditions.load_velocity_boundary_conditions!Method
load_velocity_boundary_conditions!(velocity_field, boundary_specs::Dict)

Load velocity boundary conditions from various sources.

Arguments

  • velocity_field: SHTnsVelocityField structure
  • boundary_specs: Dictionary specifying boundary sources

Examples

# No-slip boundaries at both surfaces
boundary_specs = Dict(
    :inner => (:no_slip, 0.0),
    :outer => (:no_slip, 0.0)
)

# Stress-free boundaries
boundary_specs = Dict(
    :inner => (:stress_free, 0.0),
    :outer => (:stress_free, 0.0)
)

# NetCDF file for inner, no-slip for outer
boundary_specs = Dict(
    :inner => "cmb_velocity.nc",
    :outer => (:no_slip, 0.0)
)
source
Geodynamo.BoundaryConditions.magnetic_to_qst_coefficientsMethod
magnetic_to_qst_coefficients(B_r, B_theta, B_phi, config)

Convert physical magnetic field components to QST spectral coefficients.

For magnetic fields:

  • Q: Radial component coefficients (B_r)
  • S: Spheroidal tangential component coefficients
  • T: Toroidal tangential component coefficients
source
Geodynamo.BoundaryConditions.velocity_to_qst_coefficientsMethod
velocity_to_qst_coefficients(v_r, v_theta, v_phi, config)

Convert physical velocity components (vr, vθ, v_φ) to QST spectral coefficients by using proper SHTnsKit decomposition.

The QST decomposition used by SHTnsKit:

  • Q: Radial component coefficients (transforms like scalar field)
  • S: Spheroidal horizontal component coefficients (curl-free part)
  • T: Toroidal horizontal component coefficients (divergence-free part)

Arguments

  • v_r: Radial velocity component [nlat, nlon]
  • v_theta: Colatitude velocity component [nlat, nlon]
  • v_phi: Azimuthal velocity component [nlat, nlon]
  • config: SHTnsKit configuration

Returns

  • (Q_coeffs, S_coeffs, T_coeffs): QST spectral coefficients
source

Initial Conditions

Geodynamo.InitialConditionsModule
InitialConditions

Module for loading and generating initial conditions for geodynamo simulations. Supports loading from NetCDF files, generating random fields, and setting prescribed analytical patterns.

source
Geodynamo.InitialConditions.generate_random_initial_conditions!Method
generate_random_initial_conditions!(field, field_type::Symbol;
                                   amplitude=1.0, modes_range=1:10,
                                   seed=nothing)

Generate random initial conditions for any field type.

Arguments

  • field: Field structure to initialize
  • field_type: Type of field (:temperature, :magnetic, :velocity, :composition)
  • amplitude: Overall amplitude of random perturbations
  • modes_range: Range of spherical harmonic modes to excite
  • seed: Random seed for reproducibility (optional)

Examples

# Random temperature field
generate_random_initial_conditions!(temp_field, :temperature, amplitude=0.1)

# Random magnetic field with specific modes
generate_random_initial_conditions!(mag_field, :magnetic,
                                   amplitude=0.01, modes_range=1:20, seed=42)
source
Geodynamo.InitialConditions.load_initial_conditions!Method
load_initial_conditions!(field, field_type::Symbol, file_path::String)

Load initial conditions from NetCDF file for any field type.

Arguments

  • field: Field structure (temperature, magnetic, velocity, or composition type)
  • field_type: Field type (:temperature, :magnetic, :velocity, :composition)
  • file_path: Path to NetCDF file containing initial conditions

File Format

NetCDF files should contain:

  • For scalar fields: spectral coefficients array
  • For vector fields: toroidal and poloidal spectral coefficients
  • Coordinate arrays: lm indices, radial grid
source
Geodynamo.InitialConditions.randomize_scalar_field!Method
randomize_scalar_field!(field; amplitude, lmax, domain=nothing)

Populate a scalar spectral field (temperature/composition) with random perturbations up to degree lmax. If a radial domain is provided and includes r=0, ball regularity is enforced.

source
Geodynamo.InitialConditions.save_initial_conditionsMethod
save_initial_conditions(field, field_type::Symbol, file_path::String)

Save current field state as initial conditions to NetCDF file.

This function is useful for saving generated or computed initial conditions for later use in simulations.

source
Geodynamo.InitialConditions.set_analytical_initial_conditions!Method
set_analytical_initial_conditions!(field, field_type::Symbol, pattern::Symbol;
                                  amplitude=1.0, parameters...)

Set analytical initial conditions based on predefined patterns.

Patterns

  • :conductive - Conductive temperature profile
  • :dipole - Dipolar magnetic field
  • :convective - Small convective velocity pattern
  • :stratified - Stratified composition profile

Examples

# Conductive temperature profile
set_analytical_initial_conditions!(temp_field, :temperature, :conductive)

# Earth-like dipolar magnetic field
set_analytical_initial_conditions!(mag_field, :magnetic, :dipole, amplitude=1.0)
source

Spherical Shell Geometry

Geodynamo.GeodynamoShell.create_shell_hybrid_composition_boundariesMethod
create_shell_hybrid_composition_boundaries(inner_spec, outer_spec, cfg::ShellConfig; precision=Float64)

Create composition boundaries for shell geometry. Dispatches to appropriate function based on spec types:

  • Both Tuple: fully programmatic boundaries
  • String + Tuple: hybrid (one from file, one programmatic)
  • Both String: both from files
source
Geodynamo.GeodynamoShell.create_shell_hybrid_temperature_boundariesMethod
create_shell_hybrid_temperature_boundaries(inner_spec, outer_spec, cfg::ShellConfig; precision=Float64)

Create temperature boundaries for shell geometry. Dispatches to appropriate function based on spec types:

  • Both Tuple: fully programmatic boundaries
  • String + Tuple: hybrid (one from file, one programmatic)
  • Both String: both from files
source
Geodynamo.GeodynamoShell.create_shell_radial_domainFunction
create_shell_radial_domain(nr=i_N) -> RadialDomain

Create a radial domain suitable for a spherical shell using the global parameter d_rratio for inner radius ratio. This is a thin wrapper over Geodynamo.create_radial_domain.

source

Solid Ball Geometry

Geodynamo.GeodynamoBall.ball_physical_to_spectral!Method
ball_physical_to_spectral!(phys::Geodynamo.SHTnsPhysicalField,
                           spec::Geodynamo.SHTnsSpectralField)

Wrapper for transforms in a solid sphere that enforces scalar regularity at r=0 after analysis. Use this for scalar fields (temperature, composition, etc.).

source
Geodynamo.GeodynamoBall.ball_vector_analysis!Method
ball_vector_analysis!(vec::Geodynamo.SHTnsVectorField,
                      tor::Geodynamo.SHTnsSpectralField,
                      pol::Geodynamo.SHTnsSpectralField)

Wrapper for vector analysis in a solid sphere; enforces vector regularity at r=0 after transforming to spectral toroidal/poloidal.

source
Geodynamo.GeodynamoBall.create_ball_magnetic_fieldsMethod
create_ball_magnetic_fields(T, cfg::BallConfig; nr=Geodynamo.i_N)

Create magnetic fields for a solid sphere. Since a "core" split is not used in a ball, we pass the same domain for both oc and ic to reuse the core implementation.

source
Geodynamo.GeodynamoBall.create_ball_radial_domainFunction
create_ball_radial_domain(nr=i_N) -> RadialDomain

Create a radial domain for a solid sphere (inner radius = 0). Uses a cosine-stretched grid similar to the shell for compatibility, but sets the inner radius to zero and adjusts the coordinate columns to match expectations of downstream operators.

source
Geodynamo.GeodynamoBall.enforce_ball_scalar_regularity!Method
enforce_ball_scalar_regularity!(spec::Geodynamo.SHTnsSpectralField)

Enforce scalar regularity at r=0 for solid sphere: for l>0, the scalar amplitude must vanish at r=0. Sets inner radial plane to zero for all nonzero l modes (both real and imaginary parts).

source
Geodynamo.GeodynamoBall.enforce_ball_vector_regularity!Method
enforce_ball_vector_regularity!(tor_spec::Geodynamo.SHTnsSpectralField,
                                pol_spec::Geodynamo.SHTnsSpectralField)

Enforce vector-field regularity at r=0 for solid sphere. For smooth fields, both toroidal and poloidal potentials behave like r^{l+1}, so they vanish at r=0 for all l ≥ 1. Zeros the inner radial plane for l≥1.

source

For lower-level packages used internally (SHTnsKit, PencilArrays, MPI), refer to their respective documentation.