aragog package
Submodules
aragog.core module
Core classes and functions
- class aragog.core.BoundaryConditions(_parameters: Parameters, _mesh: Mesh)
Bases:
objectBoundary conditions
- Parameters:
parameters – Parameters
mesh – Mesh
- apply_flux_boundary_conditions(state: State) None
Applies the boundary conditions to the state.
- Parameters:
state – The state to apply the boundary conditions to
- apply_flux_inner_boundary_condition(state: State) None
Applies the flux boundary condition to the state at the inner boundary.
- Parameters:
state – The state to apply the boundary conditions to
- Equivalent to CORE_BC in C code.
1: Simple core cooling 2: Prescribed heat flux 3: Prescribed temperature
- apply_flux_outer_boundary_condition(state: State) None
Applies the flux boundary condition to the state at the outer boundary.
- Parameters:
state – The state to apply the boundary conditions to
- Equivalent to SURFACE_BC in C code.
1: Grey-body atmosphere 2: Zahnle steam atmosphere 3: Couple to atmodeller 4: Prescribed heat flux 5: Prescribed temperature
- apply_temperature_boundary_conditions(temperature: ndarray[tuple[int, ...], dtype[_ScalarType_co]], temperature_basic: ndarray[tuple[int, ...], dtype[_ScalarType_co]], dTdr: ndarray[tuple[int, ...], dtype[_ScalarType_co]]) None
Conforms the temperature and dTdr at the basic nodes to temperature boundary conditions.
- Parameters:
temperature – Temperature at the staggered nodes
temperature_basic – Temperature at the basic nodes
dTdr – Temperature gradient at the basic nodes
- apply_temperature_boundary_conditions_melt(melt_fraction: ndarray[tuple[int, ...], dtype[_ScalarType_co]], melt_fraction_basic: ndarray[tuple[int, ...], dtype[_ScalarType_co]], dphidr: ndarray[tuple[int, ...], dtype[_ScalarType_co]]) None
- Conforms the melt fraction gradient dphidr at the basic nodes
to temperature boundary conditions.
- Parameters:
melt_fraction – Melt fraction at the staggered nodes
melt_fraction_basic – Melt fraction at the basic nodes
dphidr – Melt fraction gradient at the basic nodes
- class aragog.core.InitialCondition(_parameters: Parameters, _mesh: Mesh, _phases: PhaseEvaluatorCollection)
Bases:
objectInitial condition
- Parameters:
parameters – Parameters
mesh – Mesh
phases – PhaseEvaluatorCollection
- get_adiabat(pressure_basic) ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- Gets an adiabatic temperature profile by integrating
the adiatiabatic temperature gradient dTdPs from the surface. Uses the set surface temperature.
- Parameters:
nodes (Pressure field on the basic) –
- Returns:
Adiabatic temperature profile for the staggered nodes
- get_linear() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Gets a linear temperature profile
- Returns:
Linear temperature profile for the staggered nodes Only works for uniform spatial mesh.
- property temperature: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
aragog.interfaces module
Interfaces
- class aragog.interfaces.MixedPhaseEvaluatorProtocol(*args, **kwargs)
Bases:
PhaseEvaluatorProtocol,ProtocolMixed phase evaluator protocol
- liquidus() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- liquidus_gradient() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- solidus() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- solidus_gradient() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- class aragog.interfaces.PhaseEvaluatorABC
Bases:
ABCPhase evaluator ABC
- dTdPs() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
TODO: Update reference to sphinx: Solomatov (2007), Treatise on Geophysics, Eq. 3.2
- dTdrs() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- abstract delta_specific_volume() float | ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- abstract density() float | ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- abstract gravitational_acceleration() float | ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- abstract heat_capacity() float | ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- kinematic_viscosity() float | ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- abstract latent_heat() float | ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- abstract melt_fraction() float | ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- pressure: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- abstract relative_velocity() float | ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- set_pressure(pressure: ndarray[tuple[int, ...], dtype[_ScalarType_co]]) None
Sets the pressure.
- set_temperature(temperature: ndarray[tuple[int, ...], dtype[_ScalarType_co]]) None
Sets the temperature.
- temperature: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- abstract thermal_conductivity() float | ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- abstract thermal_expansivity() float | ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- update() None
Updates quantities to avoid repeat, possibly expensive, calculations.
- abstract viscosity() float | ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- class aragog.interfaces.PhaseEvaluatorProtocol(*args, **kwargs)
Bases:
ProtocolPhase evaluator protocol
- dTdPs() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- dTdrs() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- delta_specific_volume() float | ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- density() float | ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- gravitational_acceleration() float | ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- heat_capacity() float | ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- kinematic_viscosity() float | ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- latent_heat() float | ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- melt_fraction() float | ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- relative_velocity() float | ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- set_pressure(pressure: ndarray[tuple[int, ...], dtype[_ScalarType_co]]) None
- set_temperature(temperature: ndarray[tuple[int, ...], dtype[_ScalarType_co]]) None
- thermal_conductivity() float | ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- thermal_expansivity() float | ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- update() None
- viscosity() float | ndarray[tuple[int, ...], dtype[_ScalarType_co]]
aragog.mesh module
Mesh
- class aragog.mesh.AdamsWilliamsonEOS(settings: _MeshParameters, basic_radii: ndarray[tuple[int, ...], dtype[_ScalarType_co]])
Bases:
EOSAdams-Williamson equation of state (EOS).
EOS due to adiabatic self-compression from the definition of the adiabatic bulk modulus:
\[\left( \frac{{d\rho}}{{dP}} \right)_S = \frac{{\rho}}{{K_S}}\]where \(\rho\) is density, \(K_S\) the adiabatic bulk modulus, and \(S\) is entropy.
- property basic_density: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Density at basic nodes
- property basic_pressure: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Pressure at basic nodes
- get_density(pressure: float | ndarray[tuple[int, ...], dtype[_ScalarType_co]]) ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Computes density from pressure:
\[\rho(P) = \rho_s \exp(P/K_S)\]where \(\rho\) is density, \(P\) is pressure, \(\rho_s\) is surface density, and \(K_S\) is adiabatic bulk modulus.
- Parameters:
pressure – Pressure
- Returns:
Density
- get_density_from_radii(radii: float | ndarray[tuple[int, ...], dtype[_ScalarType_co]]) float | ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Computes density from radii:
\[\rho(r) = \frac{\rho_s K_S}{K_S + \rho_s g (r-r_s)}\]where \(\rho\) is density, \(r\) is radius, \(\rho_s\) is surface density, \(K_S\) is adiabatic bulk modulus, and \(r_s\) is surface radius.
- Parameters:
radii – Radii
- Returns
Density
- get_effective_density(radii) ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Computes effective density using density rho(r) integration over a spherical shell bounded by radii
- Parameters:
radii – Radii array
- Returns:
Effective Density array
- get_mass_element(radii: float | ndarray[tuple[int, ...], dtype[_ScalarType_co]]) ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Computes the mass element:
\[\frac{\delta m}{\delta r} = 4 \pi r^2 \rho\]where \(\delta m\) is the mass element, \(r\) is radius, and \(\rho\) is density.
- Parameters:
radii – Radii
- Returns:
The mass element at radii
- get_mass_within_radii(radii: float | ndarray[tuple[int, ...], dtype[_ScalarType_co]]) ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Computes mass within radii:
\[m(r) = \int 4 \pi r^2 \rho dr\]where \(m\) is mass, \(r\) is radius, and \(\rho\) is density.
The integral was evaluated using WolframAlpha.
- Parameters:
radii – Radii
- Returns:
Mass within radii
- get_mass_within_shell(radii: ndarray[tuple[int, ...], dtype[_ScalarType_co]]) ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Computes the mass within spherical shells bounded by radii.
- Parameters:
radii – Radii
- Returns:
Mass within the bounded spherical shells
- get_pressure_from_radii(radii: float | ndarray[tuple[int, ...], dtype[_ScalarType_co]]) ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Computes pressure from radii:
\[P(r) = -K_S \ln \left( 1 + \frac{\rho_s g (r-r_s)}{K_S} \right)\]where \(r\) is radius, \(K_S\) is adiabatic bulk modulus, \(P\) is pressure, \(\rho_s\) is surface density, \(g\) is gravitational acceleration, and \(r_s\) is surface radius.
- Parameters:
radii – Radii
- Returns:
Pressure
- get_pressure_gradient(pressure: float | ndarray[tuple[int, ...], dtype[_ScalarType_co]]) ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Computes the pressure gradient:
\[\frac{dP}{dr} = -g \rho\]where \(\rho\) is density, \(P\) is pressure, and \(g\) is gravitational acceleration.
- Parameters:
pressure – Pressure
- Returns:
Pressure gradient
- get_radii_from_pressure(pressure: float | ndarray[tuple[int, ...], dtype[_ScalarType_co]]) ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Computes radii from pressure:
\[P(r) = \int \frac{dP}{dr} dr = \int -g \rho_s \exp(P/K_S) dr\]And apply the boundary condition \(P=0\) at \(r=r_s\) to get:
\[r(P) = \frac{K_s \left( \exp(-P/K_S)-1 \right)}{\rho_s g} + r_s\]where \(r\) is radius, \(K_S\) is adiabatic bulk modulus, \(P\) is pressure, \(\rho_s\) is surface density, \(g\) is gravitational acceleration, and \(r_s\) is surface radius.
- Parameters:
pressure – Pressure
- Returns:
Radii
- set_staggered_pressure(staggered_radii: ndarray[tuple[int, ...], dtype[_ScalarType_co]]) None
Set staggered pressure based on staggered radii.
- property staggered_effective_density: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Effective density at staggered nodes
- property staggered_pressure: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Pressure at staggered nodes
- class aragog.mesh.EOS
Bases:
ABCGeneric EOS class
- abstract basic_density() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- abstract basic_pressure() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- abstract set_staggered_pressure(staggered_radii: ndarray[tuple[int, ...], dtype[_ScalarType_co]]) None
- abstract staggered_effective_density() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- abstract staggered_pressure() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- class aragog.mesh.FixedMesh(settings: _MeshParameters, radii: ndarray[tuple[int, ...], dtype[_ScalarType_co]], mass_radii: ndarray[tuple[int, ...], dtype[_ScalarType_co]], outer_boundary: float, inner_boundary: float)
Bases:
objectA fixed mesh
- Parameters:
settings – Mesh parameters
radii – Radii of the mesh
mass_radii – Mass coordinates of the mesh
outer_boundary – Outer boundary for computing depth below the surface
inner_boundary – Inner boundary for computing height above the base
- settings
Mesh parameters
- Type:
aragog.parser._MeshParameters
- radii
Radii of the mesh
- Type:
numpy.ndarray[tuple[int, …], numpy.dtype[numpy._typing._array_like._ScalarType_co]]
- mass_radii
Mass coordinates of the mesh
- Type:
numpy.ndarray[tuple[int, …], numpy.dtype[numpy._typing._array_like._ScalarType_co]]
- outer_boundary
Outer boundary for computing depth below the surface
- Type:
float
- inner_boundary
Inner boundary for computing height above the base
- Type:
float
- area
Surface area
- delta_mesh
Delta radii in mass coordinates
- depth
Depth below the outer boundary
- height
Height above the inner boundary
- mixing_length
Mixing length
- mixing_length_squared
Mixing length squared
- mixing_length_cubed
Mixing length cubed
- number_of_nodes
Number of nodes
- volume
Volume of the spherical shells defined between neighbouring radii
- total_volume
Total volume
- property area: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Area
- property delta_mesh: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- property depth: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- property height: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- inner_boundary: float
- mass_radii: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- property mixing_length: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- property mixing_length_cubed: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- property mixing_length_squared: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- property number_of_nodes: int
- outer_boundary: float
- radii: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- settings: _MeshParameters
- property total_volume: float
- property volume: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- class aragog.mesh.Mesh(parameters: Parameters)
Bases:
objectA staggered mesh.
The basic mesh is used for the flux calculations and the staggered mesh is used for the volume calculations.
- Parameters:
parameters – Parameters
- property basic_pressure: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- d_dr_at_basic_nodes(staggered_quantity: ndarray[tuple[int, ...], dtype[_ScalarType_co]]) ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Determines d/dr at the basic nodes of a quantity defined at the staggered nodes.
- Parameters:
staggered_quantity – A quantity defined at the staggered nodes.
- Returns:
d/dr at the basic nodes
- property dxidr: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
dxi/dr at basic nodes
- eos: EOS = Field(name=None,type=None,default=<dataclasses._MISSING_TYPE object>,default_factory=<dataclasses._MISSING_TYPE object>,init=False,repr=True,hash=None,compare=True,metadata=mappingproxy({}),kw_only=<dataclasses._MISSING_TYPE object>,_field_type=None)
- get_basic_mass_coordinates_from_spatial_coordinates(basic_coordinates: ndarray[tuple[int, ...], dtype[_ScalarType_co]]) ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Computes the basic mass coordinates from basic spatial coordinates.
- Parameters:
coordinates (Basic spatial) –
- Returns:
Basic mass coordinates
- get_constant_spacing() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Constant radius spacing across the mantle
- Returns:
Radii with constant spacing as a column vector
- get_dxidr_basic() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Computes dxidr at basic nodes.
- get_planet_density(basic_coordinates: ndarray[tuple[int, ...], dtype[_ScalarType_co]]) float
Computes the planet density.
- Parameters:
coordinates (Basic spatial) –
- Returns:
Planet effective density
- get_staggered_spatial_coordinates_from_mass_coordinates(staggered_mass_coordinates: ndarray[tuple[int, ...], dtype[_ScalarType_co]]) ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Computes the staggered spatial coordinates from staggered mass coordinates.
- Parameters:
coordinates (Staggered mass) –
- Returns:
Staggered spatial coordinates
- quantity_at_basic_nodes(staggered_quantity: ndarray[tuple[int, ...], dtype[_ScalarType_co]]) ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Determines a quantity at the basic nodes that is defined at the staggered nodes.
Uses backward and forward differences at the inner and outer radius, respectively, to obtain the quantity values of the basic nodes at the innermost and outermost nodes. When using temperature boundary conditions, values at outer boundaries will be overwritten. When using flux boundary conditions, values at outer boundaries will be used to provide estimate of individual components of heat fluxes though the total heat flux is imposed.
- Parameters:
staggered_quantity – A quantity defined at the staggered nodes
- Returns:
The quantity at the basic nodes
- quantity_at_staggered_nodes(basic_quantity: ndarray[tuple[int, ...], dtype[_ScalarType_co]]) ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Determines a quantity at the staggered nodes that is defined at the basic nodes.
Staggered nodes are always located at cell centers, whatever the mesh.
- Parameters:
basic_quantity – A quantity defined at the basic nodes
- Returns:
The quantity at the staggered nodes
- property staggered_effective_density: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- property staggered_pressure: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- volume_average(staggered_quantity: ndarray[tuple[int, ...], dtype[_ScalarType_co]]) float
- class aragog.mesh.UserDefinedEOS(settings: _MeshParameters, basic_radii: ndarray[tuple[int, ...], dtype[_ScalarType_co]])
Bases:
EOSUser defined equation of state (EOS).
Pressure field and effective density field on staggered nodes provided by the user.
- property basic_density: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Effective density at basic nodes
- property basic_pressure: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Pressure at basic nodes
- set_staggered_pressure(staggered_radii: ndarray[tuple[int, ...], dtype[_ScalarType_co]]) None
Set staggered pressure based on staggered radii.
- property staggered_effective_density: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Effective density at staggered nodes
- property staggered_pressure: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Pressure at staggered nodes
aragog.output module
Output
- class aragog.output.Output(solver: Solver)
Bases:
objectStores inputs and outputs of the models.
- property conductive_heat_flux_basic: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Conductive heat flux
- property convective_heat_flux_basic: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Convective heat flux
- property core_mass: float
Core mass computed with constant density
- property dTdr: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- property dTdrs: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- property density_basic: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Density
- property gravitational_separation_heat_flux_basic: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Gravitational separation heat flux
- property heat_capacity_basic: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Heat capacity
- property heating: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Internal heat generation at staggered nodes
- property heating_dilatation: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Internal heat generation from dilatation/compression at staggered nodes
- property heating_radio: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Internal heat generation from radioactive decay at staggered nodes
- property heating_tidal: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Internal heat generation from tidal heat dissipation at staggered nodes
- property liquidus_K_staggered: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Liquidus
- property log10_viscosity_basic: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Viscosity of the basic mesh
- property log10_viscosity_staggered: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Viscosity of the staggered mesh
- property mantle_mass: float
Mantle mass computed from the AdamsWilliamsonEOS
- property mass_radii_km_basic: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Mass radii of the basic mesh in km
- property mass_radii_km_staggered: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Mass radii of the staggered mesh in km
- property mass_staggered: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Mass of each layer on staggered mesh
- property melt_fraction_basic: float | ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Melt fraction on the basic mesh
- property melt_fraction_global: float
Volume-averaged melt fraction
- property melt_fraction_staggered: float | ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Melt fraction on the staggered mesh
- property mixing_heat_flux_basic: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Convective mixing heat flux
- plot(num_lines: int = 11, figsize: tuple = (25, 10)) None
Plots the solution with labelled lines according to time.
- Parameters:
num_lines – Number of lines to plot. Defaults to 11.
figsize – Size of the figure. Defaults to (25, 10).
- property pressure_GPa_basic: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Pressure of the basic mesh in GPa
- property pressure_GPa_staggered: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Pressure of the staggered mesh in GPa
- property radii_km_basic: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Radii of the basic mesh in km
- property radii_km_staggered: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Radii of the staggered mesh in km
- property rheological_front: float
Rheological front at the last solve iteration given user defined threshold. It is defined as a dimensionless distance with respect to the outer radius.
- property shape_basic: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Shape of the basic data
- property shape_staggered: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Shape of the staggered data
- property solidus_K_staggered: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Solidus
- property solution_top_temperature: float
Solution (last iteration) temperature at the top of the domain (planet surface)
- property super_adiabatic_temperature_gradient_basic: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Super adiabatic temperature gradient
- property temperature_K_basic: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Temperature of the basic mesh in K
- property temperature_K_staggered: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Temperature of the staggered mesh in K
- property thermal_expansivity_basic: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Thermal expansivity
- property time_range: float
- property times: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Times in years
- property total_heat_flux_basic: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Conductive heat flux
- write_at_time(file_path: str, tidx: int = -1, compress: bool = False) None
Write the state of the model at a particular time to a NetCDF4 file on the disk.
- Parameters:
file_path – Path to the output file
tidx – Index on the time axis at which to access the data
compress – Whether to compress the data
aragog.parser module
Parses the configuration file and scales and stores the parameters.
- class aragog.parser.Parameters(*, boundary_conditions: _BoundaryConditionsParameters, energy: _EnergyParameters, initial_condition: _InitialConditionParameters, mesh: _MeshParameters, phase_solid: _PhaseParameters, phase_liquid: _PhaseParameters, phase_mixed: _PhaseMixedParameters, radionuclides: list[_Radionuclide], scalings: _ScalingsParameters, solver: _SolverParameters)
Bases:
objectAssembles all the parameters.
The parameters in each section are scaled here to ensure that all the parameters are scaled (non-dimensionalised) consistently with each other.
- boundary_conditions: _BoundaryConditionsParameters
- energy: _EnergyParameters
- classmethod from_file(*filenames) Self
Parses the parameters in a configuration file(s)
- Parameters:
*filenames – Filenames of the configuration data
- initial_condition: _InitialConditionParameters
- mesh: _MeshParameters
- phase_liquid: _PhaseParameters
- phase_mixed: _PhaseMixedParameters
- phase_solid: _PhaseParameters
- static radionuclide_sections(parser: ConfigParser) list[str]
Section names relating to radionuclides
Sections relating to radionuclides must have the prefix radionuclide_
- radionuclides: list[_Radionuclide]
- scalings: _ScalingsParameters
- solver: _SolverParameters
aragog.phase module
A phase defines the equation of state (EOS) and transport properties.
- class aragog.phase.CompositePhaseEvaluator(parameters: Parameters)
Bases:
PhaseEvaluatorABCEvaluates the EOS and transport properties of a composite phase.
This combines the single phase evaluators for the liquid and solid regions with the mixed phase evaluator for the mixed phase region. This ensure that the phase properties are computed correctly for all temperatures and pressures.
- Parameters:
parameters – Parameters
- dTdPs() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
TODO: Update reference to sphinx: Solomatov (2007), Treatise on Geophysics, Eq. 3.2
- dTdrs() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- delta_specific_volume() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Difference of specific volume between melt and solid
- density() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- gravitational_acceleration() float | ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- heat_capacity() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Heat capacity
- latent_heat() float | ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Latent heat of fusion
- liquidus() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- liquidus_gradient() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- melt_fraction() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Melt fraction
- relative_velocity() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Relative velocity between melt and solid
- set_pressure(pressure: ndarray[tuple[int, ...], dtype[_ScalarType_co]]) None
Sets pressure and updates quantities that only depend on pressure
- set_temperature(temperature: ndarray[tuple[int, ...], dtype[_ScalarType_co]]) None
Sets the temperature.
- solidus() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- solidus_gradient() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- thermal_conductivity() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Thermal conductivity
- thermal_expansivity() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Thermal expansivity
- update() None
Updates quantities to avoid repeat, possibly expensive, calculations.
- viscosity() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Viscosity
- class aragog.phase.ConstantProperty(name: str, *, value: float)
Bases:
PropertyProtocolA property with a constant value
- Parameters:
name – Name of the property
value – The constant value
- name
Name of the property
- Type:
str
- value
The constant value
- Type:
float
- ndim
Number of dimensions, which is equal to zero for a constant property
- Type:
int
- eval() float
- name: str
- ndim: int = 0
- value: float
- class aragog.phase.LookupProperty1D(name: str, *, value: ndarray[tuple[int, ...], dtype[_ScalarType_co]])
Bases:
PropertyProtocolA property from a 1-D lookup
- Parameters:
name – Name of the property
value – A 2-D array, with x values in the first column and y values in the second column.
- name
Name of the property
- Type:
str
- value
A 2-D array
- Type:
npt.NDArray
- ndim
Number of dimensions, which is equal to one for a 1-D lookup
- Type:
int
- eval(pressure: ndarray[tuple[int, ...], dtype[_ScalarType_co]]) ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- gradient(pressure: ndarray[tuple[int, ...], dtype[_ScalarType_co]]) ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Computes the gradient
- name: str
- ndim: int = 1
- value: npt.NDArray
- class aragog.phase.LookupProperty2D(name: str, *, value: ndarray[tuple[int, ...], dtype[_ScalarType_co]])
Bases:
PropertyProtocolA property from a 2-D lookup
- Parameters:
name – Name of the property
value – The 2-D array
- name
Name of the property
- Type:
str
- value
The 2-D array
- Type:
npt.NDArray
- ndim
Number of dimensions, which is equal to two for a 2-D lookup
- Type:
int
- eval(temperature: ndarray[tuple[int, ...], dtype[_ScalarType_co]], pressure: ndarray[tuple[int, ...], dtype[_ScalarType_co]]) ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- name: str
- ndim: int = 2
- prepare_data_for_spline(data)
Ensure your data is on a regular grid for RectBivariateSpline
- value: npt.NDArray
- class aragog.phase.MixedPhaseEvaluator(parameters: Parameters)
Bases:
PhaseEvaluatorABCEvaluates the EOS and transport properties of a mixed phase.
This only computes quantities within the mixed phase region between the solidus and the liquidus. Computing quantities outside of this region will give incorrect results.
- Parameters:
parameters – Parameters
- settings
Mixed phase parameters
- delta_density() float | ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- delta_fusion() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- delta_specific_volume() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Difference of specific volume between melt and solid
- density() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- gravitational_acceleration() float | ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- heat_capacity() float | ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Heat capacity of the mixed phase [Equation 4, Solomatov, 2007]
- latent_heat() float
Latent heat of fusion
- liquidus() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Liquidus
- liquidus_gradient() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Liquidus gradient
- melt_fraction() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Melt fraction of the mixed phase
The melt fraction is always between zero and one.
- melt_fraction_no_clip() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Melt fraction without clipping
- porosity() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Porosity of the mixed phase, that is the volume fraction occupied by the melt
- relative_velocity() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Relative velocity between melt and solid
- set_pressure(pressure: ndarray[tuple[int, ...], dtype[_ScalarType_co]]) None
Sets pressure and updates quantities that only depend on pressure
- solidus() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Solidus
- solidus_gradient() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Solidus gradient
- thermal_conductivity() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- thermal_expansivity() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- update()
Updates quantities to avoid repeat, possibly expensive, calculations.
- viscosity() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- class aragog.phase.PhaseEvaluatorCollection(parameters: dataclasses.InitVar[Parameters])
Bases:
objectA collection of phase evaluators
Creates the phase evaluators and selects the active phase based on configuration data.
- Parameters:
parameters – Parameters
- liquid
Liquid evaluator
- solid
Solid evaluator
- mixed
Mixed evaluator
- composite
Composite evaluator
- active
The active evaluator, which is defined by configuration data
- active: PhaseEvaluatorProtocol
- composite: MixedPhaseEvaluatorProtocol
- liquid: PhaseEvaluatorProtocol
- mixed: MixedPhaseEvaluatorProtocol
- parameters: dataclasses.InitVar[Parameters]
- solid: PhaseEvaluatorProtocol
- class aragog.phase.SinglePhaseEvaluator(settings: _PhaseParameters, gravitational_acceleration: PropertyProtocol)
Bases:
PhaseEvaluatorABCContains the objects to evaluate the EOS and transport properties of a phase.
- Parameters:
settings – Phase parameters
gravitational_acceleration – PropertyProtocol
- delta_specific_volume() float | ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- density() float | ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- gravitational_acceleration() float | ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- heat_capacity() float | ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- latent_heat() float
- melt_fraction() float
- relative_velocity() float
- thermal_conductivity() float | ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- thermal_expansivity() float | ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- viscosity() float | ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- aragog.phase.setup_gravitational_acceleration(parameters: Parameters)
Sets up the gravitational acceleration property.
- Parameters:
parameters – Parameters
- Returns:
PropertyProtocol
- Return type:
gravitational_acceleration
aragog.solver module
Solver
- class aragog.solver.Evaluator(_parameters: Parameters)
Bases:
objectContains classes that evaluate quantities necessary to compute the interior evolution.
- Parameters:
_parameters – Parameters
- boundary_conditions
Boundary conditions
- initial_condition
Initial condition
- mesh
Mesh
- Type:
- phases
Evaluators for all phases
- radionuclides
Radionuclides
- boundary_conditions: BoundaryConditions
- initial_condition: InitialCondition
- phases: PhaseEvaluatorCollection
- property radionuclides: list[_Radionuclide]
- class aragog.solver.Solver(param: Parameters)
Bases:
objectSolves the interior dynamics
- Parameters:
filename – Filename of a file with configuration settings
root – Root path to the flename
- filename
Filename of a file with configuration settings
- root
Root path to the filename. Defaults to empty
- parameters
Parameters
- evaluator
Evaluator
- state
State
- dTdt(time: ndarray[tuple[int, ...], dtype[_ScalarType_co]] | float, temperature: ndarray[tuple[int, ...], dtype[_ScalarType_co]]) ndarray[tuple[int, ...], dtype[_ScalarType_co]]
dT/dt at the staggered nodes
- Parameters:
time – Time
temperature – Temperature at the staggered nodes
- Returns:
dT/dt at the staggered nodes
- classmethod from_file(filename: str | Path, root: str | Path = PosixPath('.')) Self
Parses a configuration file
- Parameters:
filename – Filename
root – Root of the filename
- Returns:
Parameters
- initialize() None
Initializes the model.
- make_tsurf_event()
Creates a temperature event function for use with an ODE solver to monitor changes in the surface temperature.The event triggers when the change exceeds the threshold, allowing the solver to stop integration.
- Returns:
terminal = True: Integration stops when the event is triggered.
direction = -1: Only triggers when the function is decreasing through zero.
- Return type:
The event has the attributes
- reset() None
This function initializes the model, while keeping the previous state of the PhaseEvaluatorCollection object. This avoids multiple loads of lookup table data when running Aragog multiple times.
- property solution: OptimizeResult
The solution.
- solve() None
- property temperature_basic: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Temperature of the basic mesh in K
- property temperature_staggered: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Temperature of the staggered mesh in K
- class aragog.solver.State(parameters: dataclasses.InitVar[Parameters], _evaluator: Evaluator)
Bases:
objectStores and updates the state at temperature and pressure.
- Parameters:
parameters – Parameters
evaluator – Evaluator
- evaluator
Evaluator
- critical_reynolds_number
Critical Reynolds number
- gravitational_separation_flux
Gravitational separation flux at the basic nodes
- heating
Total internal heat production at the staggered nodes (power per unit mass)
- heating_radio
Radiogenic heat production at the staggered nodes (power per unit mass)
- heating_tidal
Tidal heat production at the staggered nodes (power per unit mass)
- heat_flux
Heat flux at the basic nodes (power per unit area)
- inviscid_regime
True if the flow is inviscid and otherwise False, at the basic nodes
- inviscid_velocity
Inviscid velocity at the basic nodes
- is_convective
True if the flow is convecting and otherwise False, at the basic nodes
- reynolds_number
Reynolds number at the basic nodes
- super_adiabatic_temperature_gradient
Super adiabatic temperature gradient at the basic nod
- temperature_basic
Temperature at the basic nodes
- temperature_staggered
Temperature at the staggered nodes
- bottom_temperature
Temperature at the bottom basic node
- top_temperature
Temperature at the top basic node
- viscous_regime
True if the flow is viscous and otherwise False, at the basic nodes
- viscous_velocity
Viscous velocity at the basic nodes
- property bottom_temperature: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- capacitance_staggered() float | ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- conductive_heat_flux() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Conductive heat flux:
\[q_{cond} = -k \frac{\partial T}{\partial r}\]where \(k\) is thermal conductivity, \(T\) is temperature, and \(r\) is radius.
- convective_heat_flux() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Convective heat flux:
\[q_{conv} = -\rho c_p \kappa_h \left( \frac{\partial T}{\partial r} - \left( \frac{\partial T}{\partial r} \right)_S \right)\]where \(\rho\) is density, \(c_p\) is heat capacity at constant pressure, \(\kappa_h\) is eddy diffusivity, \(T\) is temperature, \(r\) is radius, and \(S\) is entropy.
- property critical_reynolds_number: float
Critical Reynolds number from Abe (1993)
- dTdr() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- dilatation_heating() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Dilatation/compression heating (power per unit mass)
- Returns:
- Dilatation/compression heating (power per unit mass) at each layer of the staggered
mesh, at a given point in time.
- dphidr() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- eddy_diffusivity() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- gravitational_separation_mass_flux() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Gravitational separation mass flux:
\[j_{grav} = \rho \phi (1 - \phi) v_{rel}\]where \(\rho\) is density, \(\phi\) is melt fraction, and \(v_{rel}\) is relative velocity.
- property heat_flux: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
The total heat flux according to the fluxes specified in the configuration.
- property heating: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
The power generation according to the heat sources specified in the configuration.
- property heating_dilatation: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
The heat source through dilation/compression.
- property heating_radio: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
The radiogenic power generation.
- property heating_tidal: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
The tidal power generation.
- property inviscid_regime: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- property inviscid_velocity: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- property is_convective: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- property mass_flux: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
The total melt mass flux according to the fluxes specified in the configuration.
- mixing_mass_flux() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Mixing mass flux:
\[j_{cm} = -\rho \kappa_h \frac{\partial \phi}{\partial r}\]where \(\rho\) is density, \(\kappa_h\) is eddy diffusivity, \(\phi\) is melt mass fraction, and \(r\) is radius.
- parameters: dataclasses.InitVar[Parameters]
- phase_basic: PhaseEvaluatorProtocol
- phase_staggered: PhaseEvaluatorProtocol
- radiogenic_heating(time: float) ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Radiogenic heating (constant with radius)
- Parameters:
time – Time
- Returns:
- Radiogenic heating (power per unit mass) at each layer of the staggered
mesh, at a given point in time.
- property reynolds_number: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- property super_adiabatic_temperature_gradient: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- property temperature_basic: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- property temperature_staggered: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- tidal_heating() ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Tidal heating at each layer of the mantle.
- Parameters:
time – Time
- Returns:
- Tidal heating (power per unit mass) at each layer of the staggered
mesh, at a given point in time.
- property top_temperature: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- update(temperature: ndarray[tuple[int, ...], dtype[_ScalarType_co]], time: float | ndarray[tuple[int, ...], dtype[_ScalarType_co]]) None
Updates the state.
The evaluation order matters because we want to minimise the number of evaluations.
- Parameters:
temperature – Temperature at the staggered nodes
pressure – Pressure at the staggered nodes
time – Time
- property viscous_regime: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
- property viscous_velocity: ndarray[tuple[int, ...], dtype[_ScalarType_co]]
aragog.utilities module
Utilities
- aragog.utilities.combine_properties(weight: Any, property1: Any, property2: Any) Any
Linear weighting of two quantities.
- Parameters:
weight – The weight to apply to property1
property1 – The value of the first property
property2 – The value of the second property
- Returns:
The combined (weighted) property
- aragog.utilities.is_file(value: Any) bool
Checks if value is a file.
- Parameters:
value – Object to be checked
- Returns:
True if the value is a file, otherwise False
- aragog.utilities.is_monotonic_increasing(some_array: ndarray[tuple[int, ...], dtype[_ScalarType_co]]) bool
Returns True if an array is monotonically increasing, otherwise returns False.
- aragog.utilities.is_number(value: Any) bool
Checks if value is a number.
- Parameters:
value – Object to be checked
- Returns:
True if the value is a number, otherwise False
- aragog.utilities.profile_decorator(func)
Decorator to profile a function
- aragog.utilities.tanh_weight(value: float | ndarray[tuple[int, ...], dtype[_ScalarType_co]], threshold: float, width: float) ndarray[tuple[int, ...], dtype[_ScalarType_co]]
Computes the tanh weight for viscosity profile and smoothing.
- Parameters:
value – Value
threshold – Threshold value
width – Width of smoothing
- Returns:
weight
Module contents
Package level variables and initialises the package logger
- aragog.aragog_file_logger(console_level=20, file_level=10, log_dir='/home/docs/checkouts/readthedocs.org/user_builds/aragog/checkouts/latest/docs') Logger
Sets up console logging and file logging according to arguments.
- Returns:
A logger
- aragog.complex_formatter() Formatter
Complex formatter for logging
- Returns:
Formatter for logging
- aragog.debug_logger() Logger
Sets up debug logging to the console.
- Returns:
A logger
- aragog.simple_formatter() Formatter
Simple formatter for logging
- Returns:
Formatter for logging