API Reference

GravityToolsNext.AbstractTempoVersionType
AbstractTempoVersion

Abstract base for different TEMPO flavors (Tempo / Tempo2). Implementations must provide:

  • tempo_data_dir(v) :: String — data directory path (TEMPO* env target)
  • tempo_cmd_path(v) :: String — resolved executable path (or bare name)
  • tempo_cmd(v) :: Cmd — ready-to-run command
  • tempo_env(v) :: Dict{String,String} — env vars to pass to the process
source
GravityToolsNext.BasicResidualStatsType
BasicResidualStats

Summary for raw (non-normalized) residuals.

Fields

  • n : number of points
  • min,max : range
  • mean : arithmetic mean
  • wmean : weighted mean (if weights given; else = mean)
  • median : sample median
  • rms : √(mean(x^2)) — uncentered
  • wrms : √(∑w x^2 / ∑w) — uncentered, weighted
  • std : √(mean((x-mean)^2)) — centered
  • wstd : √(∑w (x-wmean)^2 / ∑w) — centered, weighted
  • mad : median(|x - median(x)|) — robust spread
source
GravityToolsNext.BasicTempoOutputType
BasicTempoOutput

Basic statistical metrics from a single TEMPO2 fit iteration.

Fields:

  • rms_pre_fit_residual_us::Float64 : RMS of pre-fit residuals (µs)
  • rms_post_fit_residual_us::Float64 : RMS of post-fit residuals (µs)
  • rms_tn_post_fit_residual_us::Float64 : RMS of post-fit residuals after TN plugin (µs), if present (NaN otherwise)
  • chisqr::Float64 : Total chi-square
  • nfree::Int : Degrees of freedom
  • chisqr_red::Float64 : Reduced chi-square (chisqr / nfree)
  • pre_post::Float64 : pre-fit RMS / post-fit RMS
  • number_of_fit_parameters::Int : Number of fitted parameters
  • number_of_points_in_fit::Int : Number of TOAs in fit
  • offset_value::Float64 : Best-fit offset value
  • offset_error::Float64 : Uncertainty of the offset
  • offset_e_sqrt_n::Float64 : offset_error * sqrt(n)
source
GravityToolsNext.BasicTempoTaskType
BasicTempoTask(settings::TempoRunSettings)

Single TEMPO/TEMPO2 task with parsed iterations, optional residual statistics, and optional white-noise analysis controlled by settings.analysis.

source
GravityToolsNext.CaptureOptionsType
CaptureOptions

Controls what is captured/emitted by the TEMPO engine itself (flags / process I/O).

Fields

  • write_output::Bool : capture stdout/stderr from the engine
  • write_residuals::Bool : for Tempo2 uses -write_residuals (per internal iteration); for Tempo uses -residuals (final only)
source
GravityToolsNext.ClenshawCurtisNodesType
ClenshawCurtisNodes(n; eps=1e-6)

Chebyshev–Clenshaw–Curtis nodes on u ∈ (0,1) with optional trimming eps to avoid exact endpoints {0,1}. These u-nodes are later mapped to parameter values via prior_invcdf(prior, u).

source
GravityToolsNext.CombinedTOAEntryType
CombinedTOAEntry

Merged view of one TOA with corresponding residual entries.

Fields:

  • index::Int : sequential index after header removal (from .tim)
  • timfile_line::Int : original line number in the .tim file
  • in_fit::Bool : whether this TOA is within the requested time window
  • toa::Float64 : MJD
  • freq::Float64 : frequency (MHz)
  • signal::Float64 : residuals signal column (µs)
  • red_noise::Float64 : residual - residual_tn (µs)
  • residual::Float64 : residual with red noise (µs)
  • residual_tn::Float64 : residual with TN removed (µs)
  • uncertainty_orig::Float64: original TOA uncertainty from .tim (µs)
  • uncertainty::Float64 : transformed uncertainty from residuals file (µs)
  • backend::String : backend label
  • weight::Float64 : 1 / uncertainty^2 (or 0 if uncertainty ≤ 0 or non-finite)
source
GravityToolsNext.ContourUnitType
ContourUnit(name; min=-Inf, max=Inf, contours, from_min=true)

Refine near specified contour levels (contours) of the target metric within [min, max]. contours can be a number, a tuple, or a vector of levels.

source
GravityToolsNext.ConvergenceInfoType
ConvergenceInfo

Summary of convergence diagnostics.

  • wrms, wrms_tn, chisqr — time series with last-step deltas.
  • pre_post_final — value taken strictly from the last iteration (NaN if unavailable).
  • worst_parameterWorstParam from the last successful iteration providing a fit-parameter table; nothing otherwise.
  • thresholds — dictionary of limits used by is_converged.
  • converged — final decision under the given thresholds.
source
GravityToolsNext.ConvergenceInfoMethod
ConvergenceInfo(wrms, wrms_tn, chisqr, pre_post_final, worst_parameter, thresholds)

Builds a ConvergenceInfo and computes converged. For sequences with fewer than two points, convergence falls back to the :pre_post_final criterion.

source
GravityToolsNext.ConvergenceSeriesType
ConvergenceSeries

Holds a metric's values across iterations together with absolute/relative last-step deltas. If fewer than two points are available, final_*_delta = NaN.

source
GravityToolsNext.DiffContourUnitType
DiffContourUnit(name; min=-Inf, max=Inf, diffs, contours, from_min=true)

Refine near combinations of absolute difference thresholds (diffs) and contour levels (contours) within [min, max]. diffs/contours accept a number, tuple, or vector.

source
GravityToolsNext.DiffUnitType
DiffUnit(name; min=-Inf, max=Inf, diff, from_min=true)

Refine where the absolute change Δ(metric) exceeds diff within [min, max]. If from_min is true, scanning proceeds from the lower bound.

source
GravityToolsNext.EngineOptionsType
EngineOptions

Defines low-level options for a TEMPO run.

Fields

  • tempo_version::AbstractTempoVersion : Which TEMPO flavor to use
  • flags::String : Additional command-line flags (normalized string)
  • nits::Int : Number of internal iterations (must be ≥ 1)
  • gain::Float64 : GAIN parameter controlling convergence damping (must be > 0)

Notes

  • flags should be a single string; it is normalized (trimmed) and must not contain newlines.
  • gain is validated to be strictly positive.
source
GravityToolsNext.FitParameterType
FitParameter

Represents a single model parameter reported in the TEMPO2 fit table.

Fields:

  • name::String : parameter name (e.g., "F0", "PB", ...)
  • name_symbol::Symbol : cached symbol for fast lookup
  • pre_fit::Float64 : pre-fit value
  • post_fit::Float64 : post-fit value
  • uncertainty::Float64 : 1σ uncertainty
  • difference::Float64 : post_fit - pre_fit
  • fit_flag::Bool : whether the parameter was fitted (true) or held fixed (false)
source
GravityToolsNext.GeneralTempoResultType
GeneralTempoResult

Top-level container with everything you usually want after a run.

Fields (stored)

  • iterations : all internal iterations
  • final_index : index of the last iteration (== length(iterations))
  • last_successful_index : 0 if there were no successful iterations
  • success : quick boolean flag; true iff the run is considered successful
  • status : low-level run status (:ok | :engine_failed | :parse_failed | :files_missing | :unknown)
  • convergence : convergence summary across iterations
  • metrics : compact scalar metrics for quick ranking/comparison
  • param_estimates : map :NAME => (value, uncertainty) extracted from the final iteration
  • par_file_final : written output par-file (if available)
  • subresults : optional nested results (e.g., per-epoch, per-band, grid cells)
  • subresult_type : tag describing what subresults mean (e.g., :epoch, :band, :grid)
  • metadata : extra info (paths, timings, seeds, etc.)

Virtual properties (accessible via res.final / res.last_successful):

  • final : alias for res.iterations[res.final_index]
  • last_successful : nothing if last_successful_index == 0, otherwise res.iterations[res.last_successful_index]
  • residual_stats : alias for res.final.stats
  • white_noise : alias for res.final.white_noise_fit (white_noise_fit is also available)
source
GravityToolsNext.GridAxisType
GridAxis(name; min, max, N, rule)
GridAxis(name, values::AbstractVector)

A typed grid axis with a chosen discretization rule. Call axisvalues to obtain physical coordinates and linspace to obtain the rule's linear parameterization.

Arguments

  • name::Symbol : semantic axis name (e.g. :PBDOT, :beta0).
  • min, max : axis bounds (ignored for ExplicitRule).
  • N::Integer : number of grid points (must be ≥ 1).
  • rule : one of LinRule(), LogRule([sign]), or ExplicitRule.
  • values : alternative constructor; makes an ExplicitRule from a supplied vector of nodes.

Examples

ax1 = GridAxis(:x; min=-1.0, max=1.0, N=21, rule=LinRule())
ax2 = GridAxis(:m; min=-1e-4, max=-1e-2, N=16, rule=LogRule())  # auto-sign (−)
ax3 = GridAxis(:θ, [-1.0, -0.5, 0.0, 0.5, 1.0])                 # explicit nodes
source
GravityToolsNext.GridAxisMethod
GridAxis(name::Symbol; min, max, N, rule::AbstractGridRule=LinRule())

Construct an axis on [min, max] with N nodes, using the provided rule. Supported rules: LinRule(), LogRule([sign]), or ExplicitRule(vals) (via the positional constructor GridAxis(name, vals)).

source
GravityToolsNext.GridPriorType
GridPrior

Tabulated prior on a single parameter.

Fields

  • x : strictly increasing parameter grid (e.g., θ-values)
  • pdf : prior density on x (non-negative; integrates to 1 w.r.t. x)
  • cdf : cumulative prior on x (strictly increasing, ~[0,1])
source
GravityToolsNext.GridPriorMethod
GridPrior(x::AbstractVector{<:Real}, pdf::AbstractVector{<:Real})

Build a normalized GridPrior from tabulated x and (possibly unnormalized) pdf.

  • Sorts by x if needed;
  • Clamps negative pdf to 0;
  • Computes CDF via cumulative trapezoid;
  • Normalizes both pdf and cdf to unit mass;
  • Validates the result.
source
GravityToolsNext.InputModifiersType
InputModifiers

Describes transformations applied to the inputs of a TEMPO run.

Fields

  • override_params::Vector{TempoParameter} : parameters to upsert into the .par
  • time_start::Union{Nothing,Float64} : optional lower MJD bound (inclusive)
  • time_finish::Union{Nothing,Float64} : optional upper MJD bound (inclusive)
  • couple_f1_to_ddot::Bool : if true, auto-adjust F1 when DDOT is overridden (applied during materialization)

Notes

  • If both time_start and time_finish are provided, time_start ≤ time_finish is enforced.
source
GravityToolsNext.InternalIterationOutputType
InternalIterationOutput

Container for one internal fit iteration (or pre-iteration preamble if the first iteration did not actually run but produced logs/errors).

Fields:

  • basic::Union{BasicTempoOutput,Nothing} : basic fit metrics, or nothing if unavailable
  • fit_parameters::Union{Vector{FitParameter},Nothing} : parameter table, or nothing
  • error::TempoOutputError : error status
source
GravityToolsNext.InternalIterationResultType
InternalIterationResult

A single TEMPO iteration snapshot:

Fields

  • output : parsed TEMPO output blocks (incl. basic/error/parameters)
  • residuals : raw residual rows as read from disk (or nothing if not kept)
  • stats : residual statistics group (in-fit / in-tim), if available
  • white_noise_fit : per-backend white-noise fit result (optional)
  • metadata : extra info (paths, counts, flags, time window, etc.)
source
GravityToolsNext.LinRuleType
LinRule <: AbstractGridRule

Linear spacing rule in the physical coordinate: values(ax) is an evenly spaced LinRange(min, max, N).

source
GravityToolsNext.LocalMinimaUnitType
LocalMinimaUnit(name; min=-Inf, max=Inf, from_min=true)

Refine until a local minimum is detected along variable name within the interval [min, max].

Arguments

  • name::Symbol : the variable/axis identifier.
  • min, max : bounds of the search interval.
  • from_min : if true, sweep from the lower bound (otherwise from upper).
source
GravityToolsNext.LogRuleType
LogRule(sign)
LogRule()

Logarithmic spacing rule: uniformly spaced in log10(|x|).

Arguments

  • sign::Int : +1 or -1. Controls the sign of returned physical values. Use LogRule() (no args) to enable auto sign — inferred from the axis min at sampling time.

Notes

  • Physical values are computed as sign * 10 .^ linspace(ax), where linspace(ax) is uniform in log10(|x|).
source
GravityToolsNext.LoggingOptionsType
LoggingOptions

Lightweight logging verbosity controls for the runner.

  • verbosity — one of :silent | :warn | :info | :debug (also accepts 0..3 in the constructor)
  • with_timestamps — include timestamps in log lines
source
GravityToolsNext.MaterializedJobType
MaterializedJob

Concrete, self-contained description of a single run workspace. It keeps the original TempoRunSettings (never mutated) and records the realized filesystem layout (job_root, input/output/tmp) and resolved paths to staged inputs/outputs. Internal to TempoRun; not exported. Do not construct manually.

Fields

  • settings : the original settings used to materialize the job
  • job_root : absolute job root directory
  • input_dir : directory where input files live for this job
  • output_dir : directory where final outputs are written
  • tmp_dir : working CWD for the engine (may equal job_root)
  • run_cwd : actual CWD for the engine process
  • par_in_path : resolved path to the input par file used to create _init.par
  • tim_in_path : resolved path to the staged TIM inside the job
  • par_out_path : absolute path where the final .par is placed
  • created_paths : a list of paths created during materialization
  • notes : assorted flags (e.g., :snapshot_par, :linked_tim)
source
GravityToolsNext.PriorExecutionOptionsType
PriorExecutionOptions

Execution policy for running prior-marginalization nodes.

Fields

  • mode : node chaining policy (:independent | :chained)
  • chain_direction : traversal order when chaining (:forward | :backward)
  • chain_snapshot_par: force snapshot of per-node input par when chaining (true/false)
  • scheduler : how to schedule nodes (:serial | :distributed)
  • max_workers : cap for distributed workers (ignored for :serial; 0 means "all available")
  • workdir_layout : directory layout policy (currently only :per_node)
  • node_dir_prefix : subdirectory prefix under the base task's work dir (e.g. "nodes/node_")
  • keep_node_dirs : keep (true) or delete (false) per-node directories after run
  • on_error : error handling (:stop | :skip | :collect)

Directory naming:

  • dir_name_mode : directory name style (:index_only | :with_value)
  • index_pad : zero-padding width for node index
  • value_sig : significant digits used when embedding the parameter value into the name
source
GravityToolsNext.PriorMarginalizationSettingsType
PriorMarginalizationSettings{P,N}

Settings describing prior-based marginalization over a single parameter.

Fields

  • parameter : TEMPO parameter name (e.g., :GAMMA, :XPBDOT)
  • pin_mode : how to set the parameter per node (:force | :fixed | :fit)
  • prior : prior specification (analytic / grid / sampled)
  • nodes : rule to generate θ-nodes (mapped via prior inverse CDF, if applicable)
  • likelihood_source : metric key used to build likelihood (a key in result.metrics)
  • ref_strategy : reference choice for Δ (e.g., :prior_median | :grid_argmin | :spline_argmin | :custom_value)
  • ref_value : custom reference value (required iff ref_strategy == :custom_value)
  • representative : which θ to report at top-level (e.g., :prior_median)
  • save_node_results : whether to keep subresults (per-node GeneralTempoResult)
  • exec_options : execution policy (chaining, scheduler, dirs, error handling, naming)
  • metrics_hook : optional callback to extend the metrics dictionary

Notes

  • metrics_hook is expected to have signature like: (final_it::InternalIterationResult, conv::ConvergenceInfo, metrics::Dict{Symbol,Float64}) -> Nothing and can mutate metrics in-place to add custom scalars.
source
GravityToolsNext.PriorMarginalizedTempoTaskType
PriorMarginalizedTempoTask{T,P,N} <: SingleTempoTask

A task wrapper that evaluates a base single-tempo task over a grid of prior values θ and aggregates node results into one posterior-aware result.

Type parameters

  • T — a concrete SingleTempoTask (e.g. BasicTempoTask)
  • P — a prior specification implementing the AbstractPriorSpec API
  • N — a node rule implementing the AbstractNodeRule API

Fields

  • base_task — the underlying single task to run at each grid node
  • settings — prior definition, node grid rule, execution options, etc.
source
GravityToolsNext.RefinementSettingsType
RefinementSettings(units...; desired_refinement_level, parallel=false, params_to_save=())

Container for adaptive‑grid refinement configuration.

Arguments

  • units... : one or more refinement units (order matters).
  • desired_refinement_level : how many refinement rounds to perform (≥ 0).
  • parallel : whether the grid engine may run points in parallel.
  • params_to_save : tuple of metric names to extract at each point.
source
GravityToolsNext.RelDiffUnitType
RelDiffUnit(name; min=-Inf, max=Inf, rel_diff, from_min=true)

Refine where the relative change exceeds rel_diff within [min, max].

source
GravityToolsNext.RetentionOptionsType
RetentionOptions

Controls what data we keep in the in-memory/result structures after a run.

Fields

  • save_internal_iterations::Bool : retain all intermediate iteration results
  • save_residuals::Bool : retain residuals arrays inside results (can be large)
source
GravityToolsNext.RunArtifactsType
RunArtifacts

Filesystem view of a concrete run. All paths are absolute and can be used by task-level code to read, archive, or prune artifacts. Constructed by run_tempo_parsed.

Fields

  • job_root : job root directory
  • run_cwd : actual engine CWD (tmp/ for split layout by default)
  • input_dir : directory with staged inputs for the run
  • output_dir : directory where final outputs reside
  • tim_path : the staged TIM used by the run
  • par_out_path : final .par path (may or may not exist on failure)
  • out_path : stdout capture file path (when enabled), otherwise nothing
  • residual_paths : per-iteration residual file paths (nothing if absent)
source
GravityToolsNext.RunPathsType
RunPaths

Holds paths for a TEMPO run.

Fields

  • work_dir::String : Absolute working directory for the run
  • par_input::String : Input .par file name or relative path (resolved against work_dir; absolute paths are rejected)
  • par_output::String : Output .par file name (filename only; no directories)
  • tim_input::String : .tim file name or relative path (resolved against work_dir; absolute paths are rejected)

Notes

  • A keyword constructor RunPaths(; work_dir, par_input, tim_input, par_output=default_par_output(par_input)) is available.
source
GravityToolsNext.SampledPriorMethod
SampledPrior(path::AbstractString;
             column::Union{Int,Symbol}=1,
             comment::AbstractString="#",
             delim=nothing,
             transform::Function=identity,
             bandwidth::Union{Nothing,Real}=nothing) -> SampledPrior

Convenience constructor: read numeric samples from path and build a SampledPrior.

  • column: 1-based column index or :auto (first parseable number per line).
  • comment: lines starting with this string are skipped.
  • delim: custom delimiter for split (default: regex for commas/whitespace).
  • transform: applied to each parsed value (e.g. unit conversion).
  • bandwidth: optional KDE bandwidth override used at materialization time.
source
GravityToolsNext.TempoOutputErrorType
TempoOutputError

Describes an error encountered while parsing one iteration section.

Fields:

  • error_type::Symbol : :none if no error; otherwise a specific error tag (e.g., :missing_rms, :crash, ...)
  • message::String : diagnostic message
source
GravityToolsNext.TempoParFileType
TempoParFile

A representation of a TEMPO .par file, containing:

  • name : file name (e.g., "B1937.par")
  • dir : directory path
  • path : full path (joinpath(dir, name))
  • params : dictionary of parameters (Symbol => TempoParameter)
  • order : ordered list of parameter symbols to preserve file layout
source
GravityToolsNext.TempoParameterType

Convenient constructor. Accepts AbstractString for name (e.g. SubString), stores as String. Promotes any Real value/uncertainty to BigFloat.

source
GravityToolsNext.TempoParameterType
TempoParameter

Parameter from a TEMPO/TEMPO2 .par file.

Fields

  • name::String — canonical parameter name (may include spaces, e.g. "JUMP -be SCAMP")
  • name_symbol::Symbol — cached symbol version of name
  • value::Union{Int64, BigFloat, String, Nothing}
  • flag::Union{Int64, Nothing} — fit flag (typically -1, 0, 1), optional
  • uncertainty::Union{BigFloat, Nothing} — 1σ uncertainty, optional
source
GravityToolsNext.TempoResidualEntryType
TempoResidualEntry

One line from residuals.dat (TEMPO/TEMPO2). Values are converted to microseconds where appropriate (signal, residual, residual_tn, uncertainty).

source
GravityToolsNext.TempoRunOutputType
TempoRunOutput

Low-level result of invoking TEMPO/TEMPO2 for a single settings bundle. This is not the scientific analysis; it only captures parsed iterations, file artifacts and a coarse success/status for the engine run.

Fields

  • parsed : Vector{InternalIterationOutput} parsed from stdout
  • artifacts : RunArtifacts with absolute paths to all on-disk products
  • success : true when engine_ok ∧ parse_ok ∧ files_ok
  • status : :ok | :engine_failed | :parse_failed | :files_missing
  • n_iter : number of parsed iterations
  • exit_code : engine exit code (0 ok, 1 exception)
  • stderr_tail: last ~1000 characters of stderr, or nothing
  • started_at, finished_at, duration_s : coarse timing info
source
GravityToolsNext.TempoRunSettingsType
TempoRunSettings

Top-level settings for a single TEMPO run. Composed from smaller structures:

  • paths::RunPaths
  • engine::EngineOptions
  • modifiers::InputModifiers
  • capture::CaptureOptions
  • retention::RetentionOptions
  • analysis::WhiteNoiseOptions
  • workspace::WorkspaceOptions
  • logging::LoggingOptions
source
GravityToolsNext.TempoRunSettingsMethod
TempoRunSettings(; kwargs...) -> TempoRunSettings

Keyword-style construction with sensible defaults. Raises an error if nits < 1.

Keyword arguments

  • Paths: work_dir, par_input, par_output, tim_input
  • Execution: tempo_version, flags, nits, gain
  • Modifiers: override_params, time_start, time_finish, couple_f1_to_ddot
  • Capture: write_output, write_residuals
  • Retention: save_internal_iterations, save_residuals
  • White-noise analysis (convenience): white_noise_enabled, white_noise_scope (:final | :all) or pass analysis::WhiteNoiseOptions
  • Workspace (convenience): work_mode, job_name, overwrite (:error|:reuse|:unique|:clean), layout, temp_dir, link_tim, snapshot_par, cleanup_before_run, keep_tmp_on_success, keep_tmp_on_error, timeout_s, write_manifest, manifest_style or pass workspace::WorkspaceOptions
  • Logging (convenience): verbosity, with_timestamps or pass logging::LoggingOptions
source
GravityToolsNext.TimTOAEntryType
TimTOAEntry

Single TOA from a .tim file, after removing comments and headers.

Fields:

  • index::Int : sequential index after removing comments/headers (1-based)
  • timfile_line::Int : original line number in the .tim file (1-based)
  • toa::Float64 : MJD
  • freq::Float64 : observing frequency (MHz)
  • uncertainty::Float64: TOA uncertainty (microseconds)
  • backend::String : backend label (extracted from -be <name> if present, otherwise "unknown")
source
GravityToolsNext.WhiteNoiseBackendFitResultType
WhiteNoiseBackendFitResult

Per-backend white-noise fit result.

Fields

  • backend : backend id as Symbol
  • efac,equad : white-noise parameters used to scale uncertainties
  • offset : constant offset applied before normalization
  • ad_objective : value of the optimized objective (e.g., AD A²)
  • stats : normalized residual stats with fitted efac/equad/offset
  • success : parameters are finite/usable; diagnostics/stats computed
  • converged : optimizer reported convergence (useful even if success=true)

Notes

  • When success == false, stats contains a neutral placeholder (from _empty_norm_stats()), and the pretty-printer skips printing the stats block for this backend.
source
GravityToolsNext.WhiteNoiseFitResultType
WhiteNoiseFitResult

White-noise fit for all backends + global summary over concatenated normalized residuals.

Fields

  • by_backend : map backend => WhiteNoiseBackendFitResult
  • global_stats : summary over successful backends only
  • failed_backends : vector of backend ids that did not produce usable parameters (success=false)
source
GravityToolsNext.WhiteNoiseOptionsType
WhiteNoiseOptions(; enabled=false, scope=:final)

Options for white-noise analysis:

  • enabled — toggle analysis;
  • scope:final (last internal iteration) or :all (every internal iteration).
source
GravityToolsNext.WorkspaceOptionsType
WorkspaceOptions

Workspace/runtime controls for a single TEMPO run. The options are grouped by intent:

Lifecycle & isolation

  • work_mode::Symbol : :inplace | :jobdir
    • :inplace — run directly inside work_dir (job root = work_dir).
    • :jobdir — create/use a subdirectory of work_dir as an isolated job root.
  • job_name::Union{Nothing,String} : subdirectory name when work_mode=:jobdir.
    • If nothing, an auto name like job-YYYYmmdd-HHMMSS-<rand> is generated.
  • overwrite::Symbol : :error | :reuse | :unique | :clean
    • Behavior when the target job directory already exists (:unique appends -001, -002, ...).
    • :error — fail if exists; :reuse — use as-is; :clean — purge contents then use.
    • :overwrite is accepted for backward compatibility and treated as :clean (a warning is emitted).

Layout of files (inside the job root)

  • layout::Symbol : :flat | :split
    • :flat — files live directly in the job root.
    • :split — create input/, output/, tmp/ under the job root.
  • temp_dir::Union{Nothing,String}
    • If nothing, the execution cwd is the job root for :flat, or tmp/ for :split.
    • If non-nothing, use this subdirectory as the execution cwd (created if missing).

I/O mirroring

  • io_mirror::Union{Symbol,Int,Tuple{Symbol,Int}} = :none
    • Controls how the temp_dir (or job root) directory structure is mirrored for I/O.
    • :none (flat), :full (mirror full temp_dir), Int (first N segments), (:depth_minus, k) (first depth-k segments).
    • Ignored when layout = :flat.

Inputs capture into job root

  • link_tim::Bool : symlink the .tim file instead of copying (falls back to copy if symlink fails).
  • snapshot_par::Bool : copy the input .par for reproducibility.

Cleanup & retention

  • cleanup_before_run::Bool : remove engine artifacts in the run directory before launch.
  • keep_tmp_on_success::Bool : keep tmp/ after a successful run (useful for debugging).
  • keep_tmp_on_error::Bool : keep tmp/ if the run failed.

Runtime & manifest

  • timeout_s::Union{Nothing,Int} : kill the run after this many seconds (nothing = unlimited).
  • write_manifest::Bool : write a small manifest file in the job root.
  • manifest_style::Symbol : :json | :toml.
source
GravityToolsNext.WorstParamType
WorstParam

Named tuple describing the least stable fitted parameter: (name::Symbol, delta::Float64, uncertainty::Float64, ratio::Float64) where ratio = delta / uncertainty.

source
GravityToolsNext._PriorNodeResultType
_PriorNodeResult

Internal container holding data for a single grid node.

Fields

  • index — 1-based node position on the grid
  • theta — θ value at the node
  • workdir — absolute path to the node's working directory
  • like_metric — scalar metric used as a likelihood proxy (e.g. metrics[:chi2_fit])
  • prior_density — prior pdf value π(θ) at this node
  • result — node's GeneralTempoResult
source
Base.getindexMethod
getindex(output::InternalIterationOutput, key::Union{String,Symbol})

Convenience lookup. If key matches a field in BasicTempoOutput, returns it. Otherwise, if it matches a parameter name, returns the corresponding FitParameter. Throws KeyError if nothing matches.

source
Base.showMethod

Pretty-print Tempo/Tempo2. Supports inline one-liner via IOContext key :inline => true.

source
Base.showMethod
Base.show(io::IO, ::MIME"text/plain", grid::AdaptiveRefinement2DGrid)

Pretty-print the AdaptiveRefinement2DGrid with variable, parameter, axis, and refinement information.

source
Base.showMethod
Base.show(io::IO, par_file::TempoParFile)

Pretty-prints the .par file name, directory, and all parameters in order.

source
GravityToolsNext._apply_coupling_f1_to_ddot!Method
_apply_coupling_f1_to_ddot!(pf, overrides) -> TempoParFile

If DDOT (or F2) is present in overrides and F1 is not overridden, adjust F1 according to the heuristic ΔF1 = - F0 * (DDOTnew - DDOTbase) / D, using values from pf:

  • F0, F1 read from the current pf (pre-override state),
  • DDOT_base from pf[:DDOT] if present, otherwise pf[:F2], else 0,
  • D from pf[:D] if present, else 1.

The new F1 is upserted into pf.

source
GravityToolsNext._assemble_prior_resultMethod
_assemble_prior_result(task, thetas, nodes) -> GeneralTempoResult

Assemble the top-level prior-marginalized result:

  • choose a reference (θref, mref) and convert node metrics into log-likelihoods;
  • compute discrete posterior weights on the node grid and summarize them;
  • integrate a continuous (spline-based) posterior for moments/evidence/χ² summaries;
  • select a representative node result for the structural fields;
  • marginalize model parameter estimates across nodes using posterior weights;
  • merge metrics/extras and metadata; optionally attach node subresults.
source
GravityToolsNext._bad_prior_resultMethod
_bad_prior_result(task, bad_node, nodes_done) -> GeneralTempoResult

Return the GeneralTempoResult of a failing node augmented with prior/marginalization context in metadata. Used when on_error == :stop.

source
GravityToolsNext._check_nodes_for_errorMethod
_check_nodes_for_error(nodes, on_error) -> Union{Nothing, Int}

Return the 1-based index of the first failing node according to the policy, or nothing if execution should proceed. Currently supports on_error == :stop.

source
GravityToolsNext._empty_norm_statsMethod
_empty_norm_stats() -> NormalizedResidualStats

Return a neutral/empty NormalizedResidualStats so downstream code doesn't crash on empty groups. Internal helper.

source
GravityToolsNext._group_by_backendMethod
_group_by_backend(entries; in_fit=nothing) -> Dict{Symbol, Vector{CombinedTOAEntry}}

Group entries by backend (Symbol). Optionally filter by in_fit:

  • nothing : no filtering
  • true : only in-fit
  • false : only not in-fit

Internal helper.

source
GravityToolsNext._header_lengthMethod
_header_length(clean_lines)

Internal: returns the number of initial non-comment lines to treat as header. Backward-compatible default: first 2 non-comment lines.

source
GravityToolsNext._integrate_posteriorMethod
_integrate_posterior(thetas, metrics, prior; kwargs...) -> NamedTuple

Spline-approximate χ²(θ) from node metrics, combine with the prior to form a continuous posterior p(θ) ∝ exp(logL(θ))·π(θ), and evaluate evidence (log-domain), posterior moments, and χ² summaries.

Returns a named tuple with fields:

  • logZ (log-evidence), evidence (may be Inf if exp(logZ) overflows), mean, var, std, median, mode, quad_error,
  • θgrid, pgrid (normalized density on the returned grid),
  • chi2_post_mean_spline, chi2_at_post_mean_spline, chi2_marginalized.

Notes

  • Numerical stabilization uses c = max(Lgrid) and integrates Z0 = ∫ exp(L(θ)−c) dθ. Then: logZ = log(Z0) + c. We avoid forming Z = Z0*exp(c) directly.
  • Moments are computed as ratios in the stabilized domain: E[θ] = m1/Z0, E[θ²] = m2/Z0, Var[θ] = E[θ²] − E[θ]².
  • chi2_marginalized = −2*logZ + mref is stable even when exp(logZ) would overflow.
source
GravityToolsNext._load_final_par_fileMethod
_load_final_par_file(par_out_path::AbstractString) :: Union{TempoParFile,Nothing}

Attempt to read the final output .par file produced by TEMPO/TEMPO2. Returns a TempoParFile on success or nothing when the file is absent or cannot be parsed.

source
GravityToolsNext._make_node_taskMethod
_make_node_task(base_task, s, i, θ) -> (node_task, node_tag, node_dir, par_out)

Construct a per-node task by:

  • deriving a leaf node_dir name (no leading "nodes/"),
  • deriving node_tag = basename(node_dir),
  • deriving per-node par_output name,
  • cloning base_task with updated temp_dir, par_output, io_mirror, and an upserted parameter (s.parameter => θ) with pin-mode flag.
source
GravityToolsNext._marginalize_param_estimatesMethod
_marginalize_param_estimates(nodes, w) -> Dict{Symbol, ParamEstimate}

Compute posterior-marginalized parameter estimates across nodes using weights w. For each parameter name seen in any node, compute the pooled mean and the total uncertainty via the law of total variance.

source
GravityToolsNext._normalized_shifted_residuals_for_backendMethod
_normalized_shifted_residuals_for_backend(residuals, uncertainties_orig, efac, equad, offset)
    -> Vector{Float64}

Transform uncertainties (σ' = sqrt((efac·σ)^2 + equad^2)) and compute shifted, normalized residuals (residuals - offset) ./ σ' for a single backend. Returns a vector aligned with the inputs (no filtering performed here). Inputs are assumed to be in microseconds (µs); the result is unitless.

source
GravityToolsNext._posterior_summaryMethod
_posterior_summary(thetas, w) -> Dict{Symbol,Float64}

Return basic discrete summaries on the node grid: posterior mean, median, mode, and a normalization check.

source
GravityToolsNext._posterior_weightsMethod
_posterior_weights(thetas, loglikes, prior) -> Vector{Float64}

Compute normalized posterior weights on the node grid: wi ∝ exp(logLi) * π(θ_i), normalized with log-sum-exp stabilization.

source
GravityToolsNext._prepare_node_workdirMethod
_prepare_node_workdir(base_dir, exec, idx, param, theta) -> String

Create (if needed) and return the node working directory under base_dir. Also writes a small node_meta.txt with index/parameter/θ for traceability.

source
GravityToolsNext._resolve_referenceMethod
_resolve_reference(s, thetas, metrics) -> (θref, mref)

Choose a reference θ and its metric value used for χ² shifting and reporting. Supported strategies:

  • :grid_argmin — θ at argmin of node metrics
  • :spline_argmin — currently uses node argmin (spline argmin may be added later)
  • :prior_median — θ at prior median (nearest node)
  • :custom_value — user-provided s.ref_value
source
GravityToolsNext._rule_linvalsMethod
_rule_linvals(rule::ExplicitRule, min::Float64, max::Float64, N::Int)

For an explicit axis, the linear parameterization equals the stored values.

source
GravityToolsNext._rule_linvalsMethod
_rule_linvals(rule::LogRule, min::Float64, max::Float64, N::Int)

Return uniformly spaced samples in log10(|x|) for a logarithmic axis. Throws if min/max are zero or have different signs.

source
GravityToolsNext._rule_valuesMethod
_rule_values(rule::ExplicitRule, lin::Vector{Float64}; min_hint::Float64=NaN)

For an explicit axis, physical values equal the stored values. min_hint is accepted for a uniform API and ignored.

source
GravityToolsNext._rule_valuesMethod
_rule_values(rule::LinRule, lin::Vector{Float64}; min_hint::Float64=NaN)

For a linear axis, physical values coincide with the linear parameterization. min_hint is accepted for a uniform API and ignored.

source
GravityToolsNext._rule_valuesMethod
_rule_values(rule::LogRule, lin::Vector{Float64}; min_hint::Float64=1.0)

Map logarithmic linear samples back to physical space, using the rule sign if provided or the sign inferred from min_hint (axis min) for auto-sign.

source
GravityToolsNext._run_nodes_serial!Method
_run_nodes_serial!(task, thetas) -> Vector{_PriorNodeResult}

Run node tasks sequentially (optionally in chained mode) and collect node results. Handles per-node task derivation, execution, metric extraction, and chaining source updates. Implements on_error == :stop policy for early stopping on node failure.

source
GravityToolsNext._split_iterationsMethod
_split_iterations(output::AbstractString) -> Vector{SubString{String}}

Split the full output into iteration sections by anchor lines, BUT keep any preamble BEFORE the first anchor as part of the FIRST section (so early errors are captured). Anchors supported: lines starting with Complete fit, RMS pre-fit residual, or ss/fs.

source
GravityToolsNext._wn_maskMethod
_wn_mask(niter::Int, scope::Symbol) :: BitVector

Construct a boolean mask that selects which iterations will undergo white-noise analysis:

  • :all — all iterations are selected
  • :final — only the last iteration is selected

Unknown scopes are treated defensively as :final.

source
GravityToolsNext.ad_objective_function!Method
ad_objective_function!(residuals_norm, residuals, uncertainties_tr) -> A²

Fill residuals_norm .= residuals ./ uncertainties_tr and return Anderson–Darling A². Mutates residuals_norm as a working buffer.

source
GravityToolsNext.ad_objective_function_fit_offset!Method
ad_objective_function_fit_offset!(residuals_shifted_norm, residuals, uncertainties_tr)
    -> (A²_min, offset*)

1D minimize A² over offset, given fixed uncertainties_tr. Uses Optim.optimize. Mutates residuals_shifted_norm as a working buffer.

source
GravityToolsNext.ad_objective_function_with_offset!Method
ad_objective_function_with_offset!(residuals_shifted_norm, residuals, uncertainties_tr, offset) -> A²

Fill residuals_shifted_norm .= (residuals .- offset) ./ uncertainties_tr and return A². Mutates residuals_shifted_norm as a working buffer.

source
GravityToolsNext.axisvaluesMethod
axisvalues(ax)

Return the physical coordinates for an axis. For LogAxis, uses the rule sign (+1/−1) or infers it from ax.min when the rule is constructed with LogRule().

The result is a copy (safe to mutate by the caller).

source
GravityToolsNext.build_convergence_infoMethod
build_convergence_info(iterations; thresholds=default_convergence_thresholds())
    -> ConvergenceInfo

Constructs convergence time series from iterations that have residual statistics. pre_post_final is taken strictly from the last iteration; if unavailable there, it is NaN. The worst_parameter is taken from the last successful iteration that provides a fit-parameter table.

When only one valid point is available, the decision converged falls back to pre_post_final.

source
GravityToolsNext.build_convergence_seriesMethod
build_convergence_series(values) -> ConvergenceSeries

Constructs absolute and relative last-step deltas for a numeric sequence. Relative delta uses a stable denominator max(abs(prev), eps()).

source
GravityToolsNext.build_core_metricsMethod
build_core_metrics(final_it::InternalIterationResult, conv::ConvergenceInfo;
                   include_tim::Bool=true,
                   include_tim_even_if_same::Bool=false)
    -> Dict{Symbol,Float64}

Produce a compact, task-agnostic set of scalar metrics from the final iteration and convergence info. Missing pieces are filled with NaN.

Included keys (always):

  • :chi2fit, :chi2rfit, :wrmsfit, :wrmstn_fit
  • :prepostfinal
  • :deltawrmstn, :delta_chi2
  • :adwhitefit (global AD after white-noise fit; NaN if not available)

Additionally (only if include_tim == true and either the in-fit and in-tim sets differ or include_tim_even_if_same == true):

  • :chi2tim, :chi2rtim, :wrmstim, :wrmstn_tim
source
GravityToolsNext.build_general_tempo_resultMethod
build_general_tempo_result(
    iterations;
    par_file_final=nothing,
    subresults=GeneralTempoResult[],
    subresult_type=nothing,
    metadata=Dict(),
    metrics_hook=nothing
)

Build a GeneralTempoResult from a list of InternalIterationResults.

  • Uses the last iteration as final (index final_index = length(iterations)).
  • Computes last_successful_index (last iteration with no error and with stats); 0 if none.
  • Computes convergence across all iterations with available stats.
  • Extracts parameter estimates from the final iteration.
  • Computes metrics via build_core_metrics(final, convergence) and optionally merges a user-supplied metrics_hook(final, convergence) dictionary (values must be Real).
source
GravityToolsNext.build_internal_iteration_resultMethod
build_internal_iteration_result(output, residual_path, tim_entries; kwargs...) -> InternalIterationResult

Build a single-iteration result from parsed TEMPO output and (optionally) residuals on disk.

Behavior

  • Does not return early on TEMPO errors; instead records the error and proceeds to attempt residual loading/combination guardedly.
  • Reads residuals via a safe wrapper; if the file is missing or unreadable, residuals are nothing.
  • Combines TIM entries with residuals using a safe combiner; on length mismatch or missing inputs, combination yields nothing.
  • Optionally performs white-noise fitting (analyze_white_noise=true) on combined in-fit data when available.
  • Records per-stage flags in metadata for downstream diagnostics.

Keywords

  • save_residuals::Bool=false — keep raw residual rows inside the result
  • analyze_white_noise::Bool=false — run per-backend white-noise fit
  • time_start::Union{Nothing,Float64}=nothing — lower MJD bound for combining
  • time_finish::Union{Nothing,Float64}=nothing — upper MJD bound for combining
source
GravityToolsNext.build_node_dirnameMethod
build_node_dirname(prefix, idx; pad=3, param=:θ, theta=nothing, mode=:index_only, sig=6) -> String

Build a per-node directory name according to the selected naming mode.

  • :index_only → uses only an index with zero-padding.
  • :with_value → additionally appends __<PARAM>=<VALUE> where <VALUE> is compact-formatted.
source
GravityToolsNext.build_residual_statistics_groupMethod
build_residual_statistics_group(entries) -> ResidualStatisticsGroup

Build statistics for in-fit and in-tim sets. If there is no time window (i.e. all entries are in-fit), the same ResidualStatisticsEntry object is reused (aliased) for both fields to avoid duplicate computation/memory. You can detect this via group.in_fit === group.in_tim.

source
GravityToolsNext.build_white_noise_fitMethod
build_white_noise_fit(entries::Vector{CombinedTOAEntry}) -> WhiteNoiseFitResult

Per-backend estimation of (efac, equad, offset) on in-fit TOAs. Robust to empty groups and solver failures. Produces per-backend normalized residual stats and a global summary. Inputs are expected in microseconds (µs) for residuals and TOA uncertainties; normalized residuals are unitless.

Behavior

  • If a backend fails (throws) or yields non-finite parameters, it is marked success=false and added to failed_backends. Its stats are set to a neutral placeholder and its contribution is excluded from global_stats.
  • global_stats are computed by concatenating normalized residuals from successful backends only.
source
GravityToolsNext.calculate_point!Method
calculate_point!(p, i_ref, j_ref, grid, grid_refined, target_function, lock_obj)

Remotely calculate the value at a point in the refined grid and update both grids with locking.

source
GravityToolsNext.cell_selectorMethod
cell_selector(i_cell, j_cell, grid, ref_unit::DiffContourUnit)

Selects cells where the cell crosses a contour and the difference exceeds a threshold.

source
GravityToolsNext.cell_selectorMethod
cell_selector(i_cell, j_cell, grid, ref_unit::LocalMinimaUnit; at_corner=false)

Selects cells containing a unique local minimum within a specified range.

source
GravityToolsNext.cell_selectorMethod
cell_selector(i_cell, j_cell, grid, ref_unit::RelDiffUnit)

Selects cells where the relative difference across the cell exceeds a threshold.

source
GravityToolsNext.chisq_statsMethod
chisq_stats(res_norm; dof=length(res_norm)) -> (χ², χ²_red)

For already normalized residuals (σ'=1 scale), χ² = ∑ r².

source
GravityToolsNext.cleanup_old_tempo_filesMethod
cleanup_old_tempo_files(job::MaterializedJob) -> Nothing

Remove common TEMPO/TEMPO2 byproducts in job.job_root (and the resolved par_out_path). This is a pre-run hygiene step and never touches user inputs or previous persisted results outside the known byproduct set.

source
GravityToolsNext.cleanup_run!Method
cleanup_run!(output::TempoRunOutput, settings::TempoRunSettings; success::Bool=output.success) -> Nothing

Ergonomic overload that accepts a TempoRunOutput. Delegates to the RunArtifacts-based method. Override success if your task-level success criterion differs from the engine-level one.

source
GravityToolsNext.cleanup_run!Method
cleanup_run!(artifacts::RunArtifacts, settings::TempoRunSettings; success::Bool=true) -> Nothing

Remove the working directory (run_cwd) when it is a dedicated subdirectory (e.g., tmp/). This is a post-run hygiene step controlled by keep_tmp_on_success / keep_tmp_on_error. It never touches input_dir or output_dir and is safe to call multiple times.

source
GravityToolsNext.combine_tim_and_residualsMethod
combine_tim_and_residuals(tim_entries, residuals; time_start=nothing, time_finish=nothing)
    -> Vector{CombinedTOAEntry}

Combine .tim entries and residuals 1:1 (by index). Optional time window (MJD) marks in_fit for each entry. Lengths must match; this function throws a DimensionMismatch exception on mismatch. Use combine_tim_and_residuals_safe for a non-throwing variant.

Notes:

  • weight = 1 / uncertainty^2 if uncertainty is positive and finite; otherwise 0.0.
  • Time window is inclusive: time_start ≤ toa ≤ time_finish when both are provided.
source
GravityToolsNext.combine_tim_and_residuals_safeMethod
combine_tim_and_residuals_safe(
    tim_entries::Union{Nothing,Vector{TimTOAEntry}},
    residuals::Union{Nothing,Vector{TempoResidualEntry}};
    time_start::Union{Nothing,Float64}=nothing,
    time_finish::Union{Nothing,Float64}=nothing,
) -> Union{Nothing, Vector{CombinedTOAEntry}}

Safe (non-throwing) wrapper around combine_tim_and_residuals. Returns nothing in any of the following cases (and logs a warning):

  • tim_entries === nothing or residuals === nothing
  • either input is empty
  • length mismatch between the two inputs

Otherwise, returns the combined vector.

source
GravityToolsNext.copy_par_fileMethod
copy_par_file(par_file; new_name=nothing, suffix=nothing, new_dir=par_file.dir, deep_copy=false)

Return a (possibly shallow) copy of par_file with optional file/dir changes.

Precedence for filename:

  1. if new_name is provided -> use it,
  2. else if suffix is provided -> append _suffix before ".par",
  3. else -> keep original name.

Ensures the result ends with .par.

source
GravityToolsNext.copy_withMethod
copy_with(s::TempoRunSettings; kwargs...) -> TempoRunSettings

Create a modified copy of s with keyword overrides. You can override either whole sub-structs (analysis, workspace, logging) or individual convenience keys (e.g. white_noise_enabled, timeout_s, verbosity, ...).

Additional override-params controls:

  • override_params_clear::Bool=false: start from an empty override list.
  • override_params_delete: names (Symbol/String or a collection of them) to remove from overrides.
  • override_params_upsert::Vector{TempoParameter}=TP[]: add/replace overrides by name.

If override_params is provided, delete/upsert are applied on top of it.

source
GravityToolsNext.cpad_displayMethod
cpad_display(s, w) -> String

Center-pad string s to width w (using display width). Left pad gets the floor, right pad the remainder. Returns a new String.

source
GravityToolsNext.default_convergence_thresholdsMethod
default_convergence_thresholds() -> Dict{Symbol,Float64}

Default limits used by is_converged:

  • :abs_wrms_tn — absolute delta of TN-weighted WRMS
  • :abs_chisqr — absolute delta of χ²
  • :pre_post_final — |prepostfinal − 1|
source
GravityToolsNext.default_par_outputMethod
default_par_output(par_input::AbstractString) -> String

Generate a default output file name (no directories) by taking basename(par_input) and replacing its .par extension with *_out.par.

  • Accepts either a bare filename or any path.
  • Still enforces that the basename ends with .par (case‑insensitive).
source
GravityToolsNext.delete_params!Method
delete_params!(params, names) -> Int

Delete all parameters whose names are in names. Returns count deleted. Preserves relative order of remaining items.

source
GravityToolsNext.estimate_white_noise_ad_with_offsetMethod
estimate_white_noise_ad_with_offset(residuals, uncertainties_orig; kwargs...)
    -> (efac, equad, offset, ad_objective, converged)

Jointly estimate (efac, equad, offset) by minimizing AD A² of normalized residuals.

  • If no initial guesses are provided, builds a small grid over equad and target std(res_norm), solves efac for each pair, fits offset, and seeds Optim.optimize.
  • Mutates working buffers internally for performance (documented below).
  • converged is the boolean flag returned by Optim.converged(...) for the outer 3D solve.

Notes

  • Inputs are expected in µs for residuals and TOA uncertainties; the objective uses unitless normalized residuals.
  • On severely degenerate inputs (e.g., zero variance or infeasible seeds), the optimizer may throw; callers are expected to guard with try/catch (as done in build_white_noise_fit) and mark the backend as failed.
source
GravityToolsNext.eval_nodesMethod
eval_nodes(prior, rule::ClenshawCurtisNodes) -> θ::Vector{Float64}

Chebyshev–Clenshaw–Curtis nodes on u∈(0,1), trimmed at max(rule.eps, PRIOR_U_EPS), then mapped to parameter space via prior_invcdf.

source
GravityToolsNext.eval_nodesMethod
eval_nodes(prior, rule::ExplicitThetaNodes) -> θ::Vector{Float64}

Return a sorted, deduplicated copy of the provided θ-values (no inv-CDF mapping).

source
GravityToolsNext.eval_nodesMethod
eval_nodes(prior, rule::QuantileNodes) -> θ::Vector{Float64}

Clamp provided quantiles to (PRIORUEPS, 1-PRIORUEPS) and map via prior_invcdf.

source
GravityToolsNext.extract_tempo_parameter_from_lineMethod
extract_tempo_parameter_from_line(line) -> TempoParameter

Parse one .par line using a simple grammar:

  • name is either 1 token, or 3 tokens if the first three tokens are non-numeric strings (e.g. "JUMP -be SCAMP", "TNEF -be SCAMP").
  • Then go: value (String/Int64/BigFloat), optional flag (Int64), optional uncertainty (BigFloat).
  • If there is one trailing token after value:
    • Int64 → it's a flag
    • BigFloat → it's an uncertainty
  • If there are two trailing tokens after value and they look like Int64 then BigFloat, treat them as flag and uncertainty respectively.
source
GravityToolsNext.find_efac_for_fixed_equad_and_std!Method
find_efac_for_fixed_equad_and_std!(residuals_norm, uncertainties_tr, residuals, uncertainties_orig,
                                   equad_fixed, std_res_norm_target) -> efac

Solve for efac such that std(residuals ./ σ′(efac, equadfixed)) ≈ stdresnormtarget.

Notes:

  • Monotone in efac (for fixed equad), so we first try to bracket a root and use Bisection().
  • Mutates residuals_norm and uncertainties_tr as work buffers; no allocations.
  • Returns NaN if the target is infeasible or a robust root was not found.
source
GravityToolsNext.find_worst_parameterMethod
find_worst_parameter(fit_params) -> WorstParam | nothing

Finds the parameter with the largest absolute ratio |Δ/σ| among parameters with fit_flag = true and finite uncertainties (σ > 0). Returns a named tuple WorstParam or nothing if no suitable parameters are present.

source
GravityToolsNext.generate_par_file_nameMethod
generate_par_file_name(base_name::String, suffix::String) -> String

Generates a new .par filename by removing the existing extension (if any), appending _suffix, and adding .par.

source
GravityToolsNext.get_convergence_metricMethod
get_convergence_metric(info, key) -> Float64

Maps a threshold key to the corresponding last-step delta (or |prepostfinal−1|). Returns NaN if the metric is not available.

source
GravityToolsNext.get_inMethod
get_in(obj, path::AbstractVector{Symbol})

Follow path through nested getproperty calls: e.g. get_in(x, [:stats, :std]) is equivalent to x.stats.std.

source
GravityToolsNext.get_par_file_representationMethod
get_par_file_representation(param) -> String

Return a fixed-width .par line: name, value, optional flag, optional uncertainty. Format numerics with up to 21 significant digits; increase if needed.

source
GravityToolsNext.get_paramMethod
get_param(params, name; default=nothing) -> Union{TempoParameter,Nothing}

Return a parameter by name or default if missing.

source
GravityToolsNext.in_fit_equals_in_timMethod
in_fit_equals_in_tim(group::ResidualStatisticsGroup) -> Bool

Returns true if the group used a single aliased entry for both fields (no time window), i.e. group.in_fit === group.in_tim.

source
GravityToolsNext.is_converged_byMethod
is_converged_by(info, keys...) -> Bool

Checks convergence only for the specified keys. Metrics that are not finite are ignored for the decision.

source
GravityToolsNext.linspaceMethod
linspace(ax)

Return the linear parameterization for an axis:

  • For LinAxis, this equals the physical coordinates.
  • For LogAxis, this is uniform in log10(|x|).
  • For ExplicitAxis, returns a copy of the stored values.

The result is a copy (safe to mutate by the caller).

source
GravityToolsNext.lpad_displayMethod
lpad_display(s, w) -> String

Pad string s on the left (using display width via Unicode.textwidth) to width w. Returns a new String. If s already fits, returns String(s) unchanged.

source
GravityToolsNext.make_accessorMethod
make_accessor(spec) -> (label::String, getter::Function, key)

Create an accessor triple from a column spec. Supported forms:

  • :field — direct field access, label "field", key :field
  • "a.b" — nested path by dots, label "a.b", key "a.b"
  • (:a, :b) / Vector{Symbol} — nested path, label "a.b", key "a.b"
  • label => path_or_function — custom header; right side can be a function (row)->val or a path-like spec as above. key is the label.
  • Function — use function both as getter and as key, label "(f)".

Throws ArgumentError for unsupported types.

source
GravityToolsNext.materializeMethod
materialize(pr::SampledPrior; n_max=4096) -> GridPrior

Build (and cache) a GridPrior via KDE of the samples. Subsequent calls reuse the cache. n_max is an upper bound on the KDE output grid size (the actual size is controlled by KernelDensity.jl).

source
GravityToolsNext.materialize_jobMethod
materialize_job(s::TempoRunSettings) -> MaterializedJob

Create a reproducible run workspace according to s.workspace and stage inputs. The returned MaterializedJob stores the original settings and the realized filesystem layout and paths.

Semantics

  • work_mode:
    • :inplace — use paths.work_dir as job_root.
    • :jobdir — create or reuse a subdirectory inside paths.work_dir.
  • overwrite (when the target jobdir already exists):
    • :error — throw;
    • :reuse — use as-is (no cleanup here; see cleanup_before_run);
    • :unique — create a suffixed sibling (…-001, …-002);
    • :clean — purge contents and use the same directory.
  • layout:
    • :flatinput_dir = output_dir = job_root;
    • :split — create input/, output/, and tmp/ (or temp_dir).
  • cleanup_before_run:
    • if tmp_dir != job_root — empty only tmp_dir;
    • otherwise — remove typical TEMPO byproducts from job_root.

Staging

  • TIM: symlink when link_tim=true (fallback to copy); no-op if already in place.
  • PAR: copy snapshot when snapshot_par=true or the source is outside job_root.

Note: engine will be invoked from the chosen CWD with relative file names.

source
GravityToolsNext.param_indexMethod
param_index(params, name) -> Union{Int,Nothing}

Find index of parameter by name (compared with p.name_symbol). Returns nothing if missing.

source
GravityToolsNext.parse_fit_parametersMethod
parse_fit_parameters(section::String, ::Type{Tempo2})
    -> (Vector{FitParameter} | nothing, TempoOutputError)

Extract the parameter table between dashed lines. If dashed separators are missing, attempt a fallback search using a header-like line.

source
GravityToolsNext.parse_internal_iteration_tempo_outputMethod
parse_internal_iteration_tempo_output(section::String, ::Type{Tempo2}) -> InternalIterationOutput

Parse a single iteration section:

  1. detect known fatal errors;
  2. parse basic fit statistics;
  3. parse fit-parameter table.

If any step fails, returns an output with error set, others nothing.

source
GravityToolsNext.parse_tempo_outputMethod
parse_tempo_output(output::String, ::Type{Tempo2}) -> Vector{InternalIterationOutput}

Split full TEMPO2 output into iteration sections, keeping any preamble BEFORE the first anchor as part of the FIRST section, and parse each section.

source
GravityToolsNext.pathlike_to_symsMethod
pathlike_to_syms(x) -> Vector{Symbol}

Normalize a column/field path specification to a vector of symbols. Supported inputs:

  • Symbol
  • AbstractString with dots (e.g. "stats.std")
  • Tuple/Vector{Symbol}

Throws ArgumentError otherwise.

source
GravityToolsNext.print_aligned_tableMethod
print_aligned_table(io, d::AbstractDict; kwargs...)

Dictionary-friendly overload. Treats values as rows, and uses keys as the left-most "name" column. Keys are converted to strings. Supports the same options as the vector overload, plus:

  • name_label: header text for the name column (default: "key").
  • sort_by = :__key__ sorts by dictionary keys (name column). Use nothing to keep insertion order, or any other valid sort_by (applied to values).

All other options (cols, formats, renderers, aligns, header_align, order) are forwarded.

source
GravityToolsNext.print_aligned_tableMethod
print_aligned_table(io, rows;
    namecol, cols;
    indent=0, header=true,
    formats=Dict(), renderers=Dict(), aligns=Dict(),
    header_align=:auto,
    sort_by=nothing, order=:asc,
    names_override=nothing, name_label_override=nothing)

Flexible, allocation-friendly table printer for a vector of rows (typically structs).

Parameters

  • namecol — how to obtain the left-most "name" column from each row. One of: Symbol, dotted String (e.g. "stats.std"), tuple (:stats,:std), or a function (row)->name.
  • cols — collection of column specs (right-hand side of the table). Each element can be:
    • Symbol — field name
    • dotted String — nested path (e.g. "stats.std")
    • tuple / Vector{Symbol} — nested path
    • label => getter — custom label with either a function or a path-like spec
    • Function — getter; label becomes "(f)"
  • formats::Dict — optional number formats (Printf) per-column, keyed by the column key (as returned by make_accessor). Defaults to "%10.6f".
  • renderers::Dict — optional custom cell renderers per-column: v -> AbstractString. Takes precedence over formats.
  • aligns::Dict — optional alignment per-column: :left | :right | :center. Defaults to :right for numeric columns, :center for Bool, otherwise :left.
  • header_align — alignment for header labels:
    • :auto (default): center for right/center columns, left for text columns
    • :left, :center, or :match (same as data alignment)
  • Sorting: sort_by can be nothing (no sorting), :__name__ (sort by resolved names), any path-like spec, or a function. order is :asc or :desc.
  • Name overrides: names_override::Vector{String} can be used to directly supply the row names (useful when the name is not obtainable from the row itself). name_label_override sets the header text for the name column.

Examples

print_aligned_table(stdout, rows;
    namecol = :backend,
    cols = (
        :efac, :equad, "stats.std", "chisqr" => (:stats, :chisqr),
        :success, :converged,
    ),
    formats = Dict(
        :efac => "%10.6f", :equad => "%10.4f", "stats.std" => "%10.6f",
        "chisqr" => "%10.4f"
    ),
    renderers = Dict(
        :success => x -> (x ? "✓" : "✗"),
        :converged => x -> (x ? "✓" : "✗"),
    ),
    aligns = Dict(:success=>:center, :converged=>:center),
    header_align = :auto,
    sort_by = :ad_objective, order = :asc,
)
source
GravityToolsNext.print_white_noise_fit_report!Method
print_white_noise_fit_report!(
    io::IO,
    backend::AbstractString,
    residuals::Vector{Float64},
    uncertainties_orig::Vector{Float64},
    uncertainties_tr::Vector{Float64},
    residuals_norm::Vector{Float64},
    residuals_shifted_norm::Vector{Float64},
    efac::Float64, equad::Float64, offset::Float64, ad_objective::Float64;
    n_toas::Int = length(residuals),
    sections::AbstractSet{Symbol} = Set([
        :header, :residuals, :uncertainties_orig, :uncertainties_tr, :normalized, :shifted_normalized, :chi2, :footer
    ]),
)

Mutating workspace variant:

  • Overwrites uncertainties_tr, residuals_norm, residuals_shifted_norm.
  • Computes only what is required by sections.
source
GravityToolsNext.print_white_noise_fit_reportMethod
print_white_noise_fit_report(
    io::IO,
    backend, residuals, uncertainties_orig,
    efac, equad, offset, ad_objective;
    sections=Set([...])
)

Non-mutating convenience: allocates local work buffers and calls print_white_noise_fit_report!.

source
GravityToolsNext.prior_cdfMethod
prior_cdf(pr, x::Real) -> Float64

Prior cumulative distribution function at x. Optional but recommended for diagnostics/plots.

source
GravityToolsNext.prior_cdfMethod
prior_cdf(pr::GridPrior, x::Real) -> Float64

Evaluate the tabulated CDF with linear interpolation. Returns cdf[1] if x < x[1] and cdf[end] if x > x[end].

source
GravityToolsNext.prior_invcdfMethod
prior_invcdf(pr, u::Real) -> Float64

Inverse-CDF mapping u ∈ (0,1) to θ. Concrete prior types must implement this scalar method. The caller may clamp u into (0,1) if appropriate.

source
GravityToolsNext.prior_invcdfMethod
prior_invcdf(pr::GridPrior, u::Real) -> Float64

Invert the tabulated CDF with linear interpolation. u is clamped into [cdf[1], cdf[end]] before inversion.

source
GravityToolsNext.read_par_file!Method
read_par_file!(par_file::TempoParFile) -> TempoParFile

Reads and parses the .par file on disk.

  • Ignores blank lines and comment lines beginning with C, c, or # (after trimming).
  • Clears and fills params and order in-place.
  • If a parameter appears multiple times, the first occurrence defines the position in order, subsequent occurrences update the value in-place.
source
GravityToolsNext.read_prior_samplesMethod
read_prior_samples(path; column=1, comment="#", delim=nothing, transform=identity) -> Vector{Float64}

Read numeric samples from a text file to build an empirical prior.

  • Skips empty lines and lines starting with comment.
  • Splits each line by delim (default: regex for commas/whitespace).
  • column can be an Int (1-based) or :auto to pick the first parseable number on each line.
  • transform is applied to each parsed value (e.g., to convert units or apply a sign flip).

Returns a vector of finite Float64 values. Errors if fewer than 2 values are found.

source
GravityToolsNext.read_residual_fileMethod
read_residual_file(path::AbstractString) -> Vector{TempoResidualEntry}

Reads residuals.dat (or similar) with 5 columns: time signal residual residual_tn uncertainty. Signal, residuals, and uncertainty are converted to microseconds (×1e6).

source
GravityToolsNext.read_tim_fileMethod
read_tim_file(path::AbstractString) -> Vector{TimTOAEntry}

Reads a .tim file, skipping comment/blank lines, preserving original line numbers, and assuming the first two non-comment lines are header/meta.

Notes:

  • Expects data lines with (at least) 4 tokens where the 2nd/3rd/4th tokens are freq, MJD, uncertainty respectively (compatible with your current format).
  • Backend is taken from the -be <name> flag if present; otherwise "unknown".
  • Unparseable lines are skipped with a warning.
source
GravityToolsNext.read_tim_file_safeMethod
read_tim_file_safe(path::AbstractString) -> Union{Nothing, Vector{TimTOAEntry}}

Safe wrapper around read_tim_file. Returns nothing if the file does not exist or if any I/O/parsing exception is thrown while reading the file. Per-line parsing issues are still logged by read_tim_file as warnings.

source
GravityToolsNext.refineMethod
refine(ax::GridAxis) -> GridAxis

Refine a 1D axis by inserting midpoints to reach 2N-1 nodes.

  • LinAxis : midpoints in physical space on [min, max].
  • LogAxis : midpoints uniform in log10(|x|) (geometric means in physical space).
  • ExplicitAxis : midpoints in physical space between consecutive user nodes.

Returns a new axis; the original is unchanged.

source
GravityToolsNext.rpad_displayMethod
rpad_display(s, w) -> String

Pad string s on the right (using display width via Unicode.textwidth) to width w. Returns a new String. If s already fits, returns String(s) unchanged.

source
GravityToolsNext.run_taskMethod
run_task(task::BasicTempoTask) -> GeneralTempoResult

Run TEMPO/TEMPO2 using task.settings, parse iteration-wise outputs, attach residual/TIM-derived statistics, optionally perform white-noise fitting, and return an aggregated GeneralTempoResult.

source
GravityToolsNext.run_taskMethod
run_task(task::PriorMarginalizedTempoTask) -> GeneralTempoResult

Execute the prior-marginalized task: build θ grid, run/supply node results, and assemble the aggregated GeneralTempoResult.

If settings.exec_options.on_error == :stop, the execution stops at the first failing node (detected by _check_nodes_for_error) and returns that node's GeneralTempoResult annotated by _bad_prior_result with :prior_mode => :stopped_on_error.

source
GravityToolsNext.run_taskMethod
run_task(task::SingleTempoTask) -> GeneralTempoResult

Execute the task and return a unified result. Every concrete SingleTempoTask must provide a method.

source
GravityToolsNext.run_tempo_parsedMethod
run_tempo_parsed(settings::TempoRunSettings) :: TempoRunOutput

Materialize the workspace, run TEMPO/TEMPO2 with relative paths, parse stdout into iteration records, and assemble RunArtifacts. No post-run deletion is performed here; call cleanup_run! from task-level code after reading the artifacts.

source
GravityToolsNext.run_tempo_rawMethod
run_tempo_raw(job::MaterializedJob) -> (stdout::String, stderr::String, exit_code::Int)

Invoke TEMPO/TEMPO2 from the job's working directory (CWD) using relative paths (-f <init.par> <tim> -outpar <name.par>). Returns captured stdout/stderr and a coarse exit code: returns 0 on normal completion, 1 on Julia-side exception.

source
GravityToolsNext.run_tempo_rawMethod
run_tempo_raw(settings::TempoRunSettings)

ERROR: This method is disabled. Materialize first: job = materialize_job(settings) and then call run_tempo_raw(job). Or use run_tempo_parsed(settings).

source
GravityToolsNext.task_stage_inputs!Method
task_stage_inputs!(task::SingleTempoTask, dest_dir::AbstractString) -> Nothing

Stage/copy any task-specific input artifacts into dest_dir so that the task can run there using name-only paths. Default implementation does nothing.

Wrappers like prior-marginalization should call this before cloning a task with a new work_dir.

source
GravityToolsNext.task_with_overridesMethod
task_with_overrides(task::BasicTempoTask, overrides; work_dir) :: BasicTempoTask

Clone the task with a new work_dir and a list of TempoParameter overrides that will be upserted into the underlying TempoRunSettings.

source
GravityToolsNext.task_with_overridesMethod
task_with_overrides(task::T, overrides::AbstractVector{TempoParameter};
                    work_dir::AbstractString) -> T where {T<:SingleTempoTask}

Return a copy of task (same concrete type) with additional/overriding parameters overrides and a reassigned working directory work_dir.

Implementations must merge overrides with any existing parameter overrides present in the task's settings, giving precedence to the new ones.

source
GravityToolsNext.task_with_paramMethod
task_with_param(task::T, name::Symbol, value::Float64, flag::Int;
                work_dir::AbstractString) -> T where {T<:SingleTempoTask}

Convenience wrapper over task_with_overrides for a single parameter. By default constructs TP(String(name), value; flag=flag) and delegates. Usually there is no need to override this method.

source
GravityToolsNext.task_workdirMethod
task_workdir(task::SingleTempoTask) -> AbstractString

Return the working directory where the task runs / writes results. Used by higher-level wrappers to organize per-node subdirectories.

source
GravityToolsNext.tempo_flagMethod
tempo_flag(pin_mode::Symbol) -> Int

Map a high-level pin mode to the TEMPO flag:

  • :force-1 (force this value; do not recompute even if derived)
  • :fixed0 (set, but do not fit)
  • :fit1 (fit this parameter)
source
GravityToolsNext.transform_uncertainties!Method
transform_uncertainties!(uncertainties_tr, uncertainties_orig, efac, equad)
    -> uncertainties_tr

In-place transformation: σ' = sqrt((efac * σ)^2 + equad^2). Mutates the first argument.

source
GravityToolsNext.update_par_file!Method
update_par_file!(pf::TempoParFile; ...; prefer_overrides=true) -> TempoParFile

Apply overrides and optionally set NITS/GAIN/START/FINISH with a clear precedence rule.

  • couple_f1_to_ddot::Bool=false — if true and DDOT is overridden (while F1 is not), auto-adjust F1 using ΔF1 ≈ -F0·ΔDDOT/D.
source
GravityToolsNext.update_points_to_calculate!Method
update_points_to_calculate!(points_to_calculate, i_ref_init, j_ref_init, grid, grid_refined, refined_points_status)

Update the list of points to calculate for a given refined point (parallel version).

source
GravityToolsNext.upsert_param!Method
upsert_param!(params, p) -> Vector{TempoParameter}

Update by name or append if not present. Keeps order. Mutates params.

Rules:

  • If parameter with this name is new → push p.
  • If exists → replace value/flag intelligently:
    • If p.value === nothing, keep the old value.
    • If p.flag === nothing, keep the old flag.
    • Otherwise take the provided values.
source
GravityToolsNext.upsert_params!Method
upsert_params!(params, new_params) -> Vector{TempoParameter}

Apply upsert_param! for each element of new_params. Order of new elements is preserved.

source
GravityToolsNext.validateMethod
validate(v::AbstractTempoVersion) -> Bool

Check that data directory exists and the executable is available (or at least named). Returns true if OK, otherwise throws.

source
GravityToolsNext.validateMethod
validate(s::TempoRunSettings) -> Bool

Alias for [validate_inputs_exist] — quick pre-materialization check that input .par/.tim exist relative to work_dir and that par_output is a bare file name. Does not create directories.

source
GravityToolsNext.validate_inputs_existMethod
validate_inputs_exist(s::TempoRunSettings) -> Bool

Validate inputs before materialization:

  • checks that par_input and tim_input exist relative to work_dir;
  • ensures par_output is a file name (no directories).

Does not create directories. Returns true if valid, otherwise throws.

source
GravityToolsNext.validate_priorMethod
validate_prior(pr::GridPrior; atol_mass=1e-3) -> Bool

Validate monotonicity and normalization of a GridPrior. Throws an error with a diagnostic message if invalid.

source
GravityToolsNext.weighted_meanMethod
weighted_mean(res) -> (μ, se_μ)
weighted_mean(res, unc) -> (μ_w, se_μ)

Unweighted fast path avoids allocating ones(N). Weighted uses w = 1/unc^2. Standard errors follow the same conventions as in the legacy printout.

source
GravityToolsNext.weighted_mean_invvarMethod
weighted_mean_invvar(residuals, uncertainties_tr) -> (μ_w, se_μ)

Inverse-variance weighted mean with w = 1/σ′^2 without allocations. Skips entries with non-finite residuals/σ′ or σ′ ≤ 0. Returns (NaN, NaN) if no valid points.

source
GravityToolsNext.weighted_stdMethod
weighted_std(res) / weighted_std(res, unc) -> (σ, se_σ)

Centered at (weighted) mean. Unweighted uses population σ (divide by N), to match the diagnostic behavior you had.

source
GravityToolsNext.write_node_metadataMethod
write_node_metadata(path; meta...) -> Nothing

Write a simple key: value metadata file at path (e.g., "node_meta.txt").

  • Real values are formatted via format_short.
  • AbstractVector{<:Real} values are comma-joined with each element formatted via format_short.
  • Other vectors are comma-joined via string.
  • Other scalars are stringified.
source
GravityToolsNext.write_par_file!Method
write_par_file!(par_file::TempoParFile) -> TempoParFile

Writes parameters to par_file.path, preserving the order defined in order. If par_file.dir does not exist, it will be created.

source